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 PetscMPIInt rank; 13c9493c35SStefano Zampini PetscBool iascii,isdraw,isglvis,isbinary; 1447c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 159a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 169a42bb27SBarry Smith PetscBool ismatlab; 179a42bb27SBarry Smith #endif 1847c6ae99SBarry Smith 1947c6ae99SBarry Smith PetscFunctionBegin; 209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)da),&rank)); 2147c6ae99SBarry Smith 229566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 239566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw)); 249566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERGLVIS,&isglvis)); 259566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary)); 269a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 279566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERMATLAB,&ismatlab)); 289a42bb27SBarry Smith #endif 2947c6ae99SBarry Smith if (iascii) { 3047c6ae99SBarry Smith PetscViewerFormat format; 3147c6ae99SBarry Smith 329566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 339566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format)); 3476a8abe0SBarry Smith if (format == PETSC_VIEWER_LOAD_BALANCE) { 3576a8abe0SBarry Smith PetscInt i,nmax = 0,nmin = PETSC_MAX_INT,navg = 0,*nz,nzlocal; 3676a8abe0SBarry Smith DMDALocalInfo info; 3776a8abe0SBarry Smith PetscMPIInt size; 389566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)da),&size)); 399566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da,&info)); 4076a8abe0SBarry Smith nzlocal = info.xm*info.ym*info.zm; 419566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size,&nz)); 429566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&nzlocal,1,MPIU_INT,nz,1,MPIU_INT,PetscObjectComm((PetscObject)da))); 4376a8abe0SBarry Smith for (i=0; i<(PetscInt)size; i++) { 4476a8abe0SBarry Smith nmax = PetscMax(nmax,nz[i]); 4576a8abe0SBarry Smith nmin = PetscMin(nmin,nz[i]); 4676a8abe0SBarry Smith navg += nz[i]; 4776a8abe0SBarry Smith } 489566063dSJacob Faibussowitsch PetscCall(PetscFree(nz)); 4976a8abe0SBarry Smith navg = navg/size; 509566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Load Balance - Grid Points: Min %D avg %D max %D\n",nmin,navg,nmax)); 5176a8abe0SBarry Smith PetscFunctionReturn(0); 5276a8abe0SBarry Smith } 538ec8862eSJed Brown if (format != PETSC_VIEWER_ASCII_VTK_DEPRECATED && format != PETSC_VIEWER_ASCII_VTK_CELL_DEPRECATED && format != PETSC_VIEWER_ASCII_GLVIS) { 54aa219208SBarry Smith DMDALocalInfo info; 559566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da,&info)); 569566063dSJacob Faibussowitsch PetscCall(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)); 579566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"X range of indices: %D %D, Y range of indices: %D %D, Z range of indices: %D %D\n", 585f80ce2aSJacob Faibussowitsch info.xs,info.xs+info.xm,info.ys,info.ys+info.ym,info.zs,info.zs+info.zm)); 5947c6ae99SBarry Smith #if !defined(PETSC_USE_COMPLEX) 606636e97aSMatthew G Knepley if (da->coordinates) { 6147c6ae99SBarry Smith PetscInt last; 6247c6ae99SBarry Smith const PetscReal *coors; 639566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(da->coordinates,&coors)); 649566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(da->coordinates,&last)); 6547c6ae99SBarry Smith last = last - 3; 669566063dSJacob Faibussowitsch PetscCall(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])); 679566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(da->coordinates,&coors)); 6847c6ae99SBarry Smith } 6947c6ae99SBarry Smith #endif 709566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 719566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 728135c375SStefano Zampini } else if (format == PETSC_VIEWER_ASCII_GLVIS) { 739566063dSJacob Faibussowitsch PetscCall(DMView_DA_GLVis(da,viewer)); 74616157d6SJed Brown } else { 759566063dSJacob Faibussowitsch PetscCall(DMView_DA_VTK(da,viewer)); 7647c6ae99SBarry Smith } 7747c6ae99SBarry Smith } else if (isdraw) { 7847c6ae99SBarry Smith PetscDraw draw; 7947c6ae99SBarry Smith PetscReal ymin = -1.0,ymax = (PetscReal)dd->N; 8047c6ae99SBarry Smith PetscReal xmin = -1.0,xmax = (PetscReal)((dd->M+2)*dd->P),x,y,ycoord,xcoord; 818ea3bf28SBarry Smith PetscInt k,plane,base; 828ea3bf28SBarry Smith const PetscInt *idx; 8347c6ae99SBarry Smith char node[10]; 8447c6ae99SBarry Smith PetscBool isnull; 855f80ce2aSJacob Faibussowitsch PetscErrorCode ierr; 8647c6ae99SBarry Smith 879566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(viewer,0,&draw)); 889566063dSJacob Faibussowitsch PetscCall(PetscDrawIsNull(draw,&isnull)); 895b399a63SLisandro Dalcin if (isnull) PetscFunctionReturn(0); 9047c6ae99SBarry Smith 919566063dSJacob Faibussowitsch PetscCall(PetscDrawCheckResizedWindow(draw)); 929566063dSJacob Faibussowitsch PetscCall(PetscDrawClear(draw)); 939566063dSJacob Faibussowitsch PetscCall(PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax)); 945b399a63SLisandro Dalcin 959566063dSJacob Faibussowitsch ierr = PetscDrawCollectiveBegin(draw);PetscCall(ierr); 9647c6ae99SBarry Smith /* first processor draw all node lines */ 97dd400576SPatrick Sanan if (rank == 0) { 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++) { 1019566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK)); 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++) { 1059566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK)); 10647c6ae99SBarry Smith } 10747c6ae99SBarry Smith } 10847c6ae99SBarry Smith } 1099566063dSJacob Faibussowitsch ierr = PetscDrawCollectiveEnd(draw);PetscCall(ierr); 1109566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 1119566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw)); 11247c6ae99SBarry Smith 1139566063dSJacob Faibussowitsch ierr = PetscDrawCollectiveBegin(draw);PetscCall(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 1239566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED)); 1249566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED)); 1259566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED)); 1269566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED)); 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 */ 1329566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(node,sizeof(node),"%d",(int)rank)); 1339566063dSJacob Faibussowitsch PetscCall(PetscDrawString(draw,xmin+(dd->M+1)*k+.2,ymin+.3,PETSC_DRAW_RED,node)); 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++) { 1389566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(node,sizeof(node),"%d",(int)base++)); 1399566063dSJacob Faibussowitsch PetscCall(PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node)); 14047c6ae99SBarry Smith } 14147c6ae99SBarry Smith } 14247c6ae99SBarry Smith 14347c6ae99SBarry Smith } 14447c6ae99SBarry Smith } 1459566063dSJacob Faibussowitsch ierr = PetscDrawCollectiveEnd(draw);PetscCall(ierr); 1469566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 1479566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw)); 14847c6ae99SBarry Smith 1499566063dSJacob Faibussowitsch ierr = PetscDrawCollectiveBegin(draw);PetscCall(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; 1559566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(da->ltogmap,&idx)); 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); 1739566063dSJacob Faibussowitsch PetscCall(PetscDrawString(draw,xcoord/dd->w,ycoord,PETSC_DRAW_BLUE,node)); 174565245c5SBarry Smith base++; 17547c6ae99SBarry Smith } 17647c6ae99SBarry Smith } 1779566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(da->ltogmap,&idx)); 17847c6ae99SBarry Smith } 17947c6ae99SBarry Smith } 1809566063dSJacob Faibussowitsch ierr = PetscDrawCollectiveEnd(draw);PetscCall(ierr); 1819566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 1829566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw)); 1839566063dSJacob Faibussowitsch PetscCall(PetscDrawSave(draw)); 184c9493c35SStefano Zampini } else if (isglvis) { 1859566063dSJacob Faibussowitsch PetscCall(DMView_DA_GLVis(da,viewer)); 1869a42bb27SBarry Smith } else if (isbinary) { 1879566063dSJacob Faibussowitsch PetscCall(DMView_DA_Binary(da,viewer)); 1889a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 1899a42bb27SBarry Smith } else if (ismatlab) { 1909566063dSJacob Faibussowitsch PetscCall(DMView_DA_Matlab(da,viewer)); 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 23047c6ae99SBarry Smith PetscFunctionBegin; 2317a8be351SBarry Smith PetscCheck(stencil_type != DMDA_STENCIL_BOX || (bx != DM_BOUNDARY_MIRROR && by != DM_BOUNDARY_MIRROR && bz != DM_BOUNDARY_MIRROR),PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Mirror boundary and box stencil"); 2329566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject) da, &comm)); 2333855c12bSBarry Smith #if !defined(PETSC_USE_64BIT_INDICES) 2347a8be351SBarry Smith PetscCheck(((PetscInt64) M)*((PetscInt64) N)*((PetscInt64) P)*((PetscInt64) dof) <= (PetscInt64) PETSC_MPI_INT_MAX,comm,PETSC_ERR_INT_OVERFLOW,"Mesh of %D by %D by %D by %D (dof) is too large for 32 bit indices",M,N,P,dof); 2353855c12bSBarry Smith #endif 2363855c12bSBarry Smith 2379566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 2389566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm,&rank)); 23947c6ae99SBarry Smith 24047c6ae99SBarry Smith if (m != PETSC_DECIDE) { 241*08401ef6SPierre Jolivet PetscCheck(m >= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m); 242*08401ef6SPierre Jolivet else PetscCheck(m <= size,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size); 24347c6ae99SBarry Smith } 24447c6ae99SBarry Smith if (n != PETSC_DECIDE) { 245*08401ef6SPierre Jolivet PetscCheck(n >= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n); 246*08401ef6SPierre Jolivet else PetscCheck(n <= size,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size); 24747c6ae99SBarry Smith } 24847c6ae99SBarry Smith if (p != PETSC_DECIDE) { 249*08401ef6SPierre Jolivet PetscCheck(p >= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Z direction: %D",p); 250*08401ef6SPierre Jolivet else PetscCheck(p <= size,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Z direction: %D %d",p,size); 25147c6ae99SBarry Smith } 2522c71b3e2SJacob Faibussowitsch PetscCheckFalse((m > 0) && (n > 0) && (p > 0) && (m*n*p != size),PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"m %D * n %D * p %D != size %d",m,n,p,size); 25347c6ae99SBarry Smith 25447c6ae99SBarry Smith /* Partition the array among the processors */ 25547c6ae99SBarry Smith if (m == PETSC_DECIDE && n != PETSC_DECIDE && p != PETSC_DECIDE) { 25647c6ae99SBarry Smith m = size/(n*p); 25747c6ae99SBarry Smith } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) { 25847c6ae99SBarry Smith n = size/(m*p); 25947c6ae99SBarry Smith } else if (m != PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) { 26047c6ae99SBarry Smith p = size/(m*n); 26147c6ae99SBarry Smith } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) { 26247c6ae99SBarry Smith /* try for squarish distribution */ 263369cc0aeSBarry Smith m = (int)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)N*p))); 26447c6ae99SBarry Smith if (!m) m = 1; 26547c6ae99SBarry Smith while (m > 0) { 26647c6ae99SBarry Smith n = size/(m*p); 26747c6ae99SBarry Smith if (m*n*p == size) break; 26847c6ae99SBarry Smith m--; 26947c6ae99SBarry Smith } 2707a8be351SBarry Smith PetscCheck(m,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad p value: p = %D",p); 27147c6ae99SBarry Smith if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;} 27247c6ae99SBarry Smith } else if (m == PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) { 27347c6ae99SBarry Smith /* try for squarish distribution */ 274369cc0aeSBarry Smith m = (int)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)P*n))); 27547c6ae99SBarry Smith if (!m) m = 1; 27647c6ae99SBarry Smith while (m > 0) { 27747c6ae99SBarry Smith p = size/(m*n); 27847c6ae99SBarry Smith if (m*n*p == size) break; 27947c6ae99SBarry Smith m--; 28047c6ae99SBarry Smith } 2817a8be351SBarry Smith PetscCheck(m,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad n value: n = %D",n); 28247c6ae99SBarry Smith if (M > P && m < p) {PetscInt _m = m; m = p; p = _m;} 28347c6ae99SBarry Smith } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) { 28447c6ae99SBarry Smith /* try for squarish distribution */ 285369cc0aeSBarry Smith n = (int)(0.5 + PetscSqrtReal(((PetscReal)N)*((PetscReal)size)/((PetscReal)P*m))); 28647c6ae99SBarry Smith if (!n) n = 1; 28747c6ae99SBarry Smith while (n > 0) { 28847c6ae99SBarry Smith p = size/(m*n); 28947c6ae99SBarry Smith if (m*n*p == size) break; 29047c6ae99SBarry Smith n--; 29147c6ae99SBarry Smith } 2927a8be351SBarry Smith PetscCheck(n,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad m value: m = %D",n); 29347c6ae99SBarry Smith if (N > P && n < p) {PetscInt _n = n; n = p; p = _n;} 29447c6ae99SBarry Smith } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) { 29547c6ae99SBarry Smith /* try for squarish distribution */ 296369cc0aeSBarry Smith n = (PetscInt)(0.5 + PetscPowReal(((PetscReal)N*N)*((PetscReal)size)/((PetscReal)P*M),(PetscReal)(1./3.))); 29747c6ae99SBarry Smith if (!n) n = 1; 29847c6ae99SBarry Smith while (n > 0) { 29947c6ae99SBarry Smith pm = size/n; 30047c6ae99SBarry Smith if (n*pm == size) break; 30147c6ae99SBarry Smith n--; 30247c6ae99SBarry Smith } 30347c6ae99SBarry Smith if (!n) n = 1; 304369cc0aeSBarry Smith m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)P*n))); 30547c6ae99SBarry Smith if (!m) m = 1; 30647c6ae99SBarry Smith while (m > 0) { 30747c6ae99SBarry Smith p = size/(m*n); 30847c6ae99SBarry Smith if (m*n*p == size) break; 30947c6ae99SBarry Smith m--; 31047c6ae99SBarry Smith } 31147c6ae99SBarry Smith if (M > P && m < p) {PetscInt _m = m; m = p; p = _m;} 312*08401ef6SPierre Jolivet } else PetscCheck(m*n*p == size,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition"); 31347c6ae99SBarry Smith 314*08401ef6SPierre Jolivet PetscCheck(m*n*p == size,PetscObjectComm((PetscObject)da),PETSC_ERR_PLIB,"Could not find good partition"); 315*08401ef6SPierre Jolivet PetscCheck(M >= m,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m); 316*08401ef6SPierre Jolivet PetscCheck(N >= n,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n); 317*08401ef6SPierre Jolivet PetscCheck(P >= p,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Partition in z direction is too fine! %D %D",P,p); 31847c6ae99SBarry Smith 31947c6ae99SBarry Smith /* 32047c6ae99SBarry Smith Determine locally owned region 32147c6ae99SBarry Smith [x, y, or z]s is the first local node number, [x, y, z] is the number of local nodes 32247c6ae99SBarry Smith */ 32347c6ae99SBarry Smith 32447c6ae99SBarry Smith if (!lx) { 3259566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m, &dd->lx)); 32647c6ae99SBarry Smith lx = dd->lx; 3278865f1eaSKarl Rupp for (i=0; i<m; i++) lx[i] = M/m + ((M % m) > (i % m)); 32847c6ae99SBarry Smith } 32947c6ae99SBarry Smith x = lx[rank % m]; 33047c6ae99SBarry Smith xs = 0; 3318865f1eaSKarl Rupp for (i=0; i<(rank%m); i++) xs += lx[i]; 3322c71b3e2SJacob Faibussowitsch PetscCheckFalse((x < s) && ((m > 1) || (bx == DM_BOUNDARY_PERIODIC)),PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local x-width of domain x %D is smaller than stencil width s %D",x,s); 33347c6ae99SBarry Smith 33447c6ae99SBarry Smith if (!ly) { 3359566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &dd->ly)); 33647c6ae99SBarry Smith ly = dd->ly; 3378865f1eaSKarl Rupp for (i=0; i<n; i++) ly[i] = N/n + ((N % n) > (i % n)); 33847c6ae99SBarry Smith } 33947c6ae99SBarry Smith y = ly[(rank % (m*n))/m]; 3402c71b3e2SJacob Faibussowitsch PetscCheckFalse((y < s) && ((n > 1) || (by == DM_BOUNDARY_PERIODIC)),PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local y-width of domain y %D is smaller than stencil width s %D",y,s); 34130729d88SBarry Smith 34247c6ae99SBarry Smith ys = 0; 3438865f1eaSKarl Rupp for (i=0; i<(rank % (m*n))/m; i++) ys += ly[i]; 34447c6ae99SBarry Smith 34547c6ae99SBarry Smith if (!lz) { 3469566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(p, &dd->lz)); 34747c6ae99SBarry Smith lz = dd->lz; 3488865f1eaSKarl Rupp for (i=0; i<p; i++) lz[i] = P/p + ((P % p) > (i % p)); 34947c6ae99SBarry Smith } 35047c6ae99SBarry Smith z = lz[rank/(m*n)]; 351bcea557cSEthan Coon 352fdc81ce8SEthan Coon /* note this is different than x- and y-, as we will handle as an important special 353bff4a2f0SMatthew G. Knepley case when p=P=1 and DM_BOUNDARY_PERIODIC and s > z. This is to deal with 2D problems 354fdc81ce8SEthan Coon in a 3D code. Additional code for this case is noted with "2d case" comments */ 3556f951b95Secoon twod = PETSC_FALSE; 3568865f1eaSKarl Rupp if (P == 1) twod = PETSC_TRUE; 3572c71b3e2SJacob Faibussowitsch else PetscCheckFalse((z < s) && ((p > 1) || (bz == DM_BOUNDARY_PERIODIC)),PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local z-width of domain z %D is smaller than stencil width s %D",z,s); 35847c6ae99SBarry Smith zs = 0; 3598865f1eaSKarl Rupp for (i=0; i<(rank/(m*n)); i++) zs += lz[i]; 36047c6ae99SBarry Smith ye = ys + y; 36147c6ae99SBarry Smith xe = xs + x; 36247c6ae99SBarry Smith ze = zs + z; 36347c6ae99SBarry Smith 36488661749SPeter Brune /* determine ghost region (Xs) and region scattered into (IXs) */ 365d9c9ebe5SPeter Brune if (xs-s > 0) { 366d9c9ebe5SPeter Brune Xs = xs - s; IXs = xs - s; 36788661749SPeter Brune } else { 3688865f1eaSKarl Rupp if (bx) Xs = xs - s; 3698865f1eaSKarl Rupp else Xs = 0; 37088661749SPeter Brune IXs = 0; 37188661749SPeter Brune } 372d9c9ebe5SPeter Brune if (xe+s <= M) { 373d9c9ebe5SPeter Brune Xe = xe + s; IXe = xe + s; 37488661749SPeter Brune } else { 37588661749SPeter Brune if (bx) { 376d9c9ebe5SPeter Brune Xs = xs - s; Xe = xe + s; 3778865f1eaSKarl Rupp } else Xe = M; 37888661749SPeter Brune IXe = M; 37988661749SPeter Brune } 38047c6ae99SBarry Smith 381bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC || bx == DM_BOUNDARY_MIRROR) { 382d9c9ebe5SPeter Brune IXs = xs - s; 383d9c9ebe5SPeter Brune IXe = xe + s; 384d9c9ebe5SPeter Brune Xs = xs - s; 385d9c9ebe5SPeter Brune Xe = xe + s; 38688661749SPeter Brune } 38788661749SPeter Brune 388d9c9ebe5SPeter Brune if (ys-s > 0) { 389d9c9ebe5SPeter Brune Ys = ys - s; IYs = ys - s; 39088661749SPeter Brune } else { 3918865f1eaSKarl Rupp if (by) Ys = ys - s; 3928865f1eaSKarl Rupp else Ys = 0; 39388661749SPeter Brune IYs = 0; 39488661749SPeter Brune } 395d9c9ebe5SPeter Brune if (ye+s <= N) { 396d9c9ebe5SPeter Brune Ye = ye + s; IYe = ye + s; 39788661749SPeter Brune } else { 3988865f1eaSKarl Rupp if (by) Ye = ye + s; 3998865f1eaSKarl Rupp else Ye = N; 40088661749SPeter Brune IYe = N; 40188661749SPeter Brune } 40288661749SPeter Brune 403bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_PERIODIC || by == DM_BOUNDARY_MIRROR) { 404d9c9ebe5SPeter Brune IYs = ys - s; 405d9c9ebe5SPeter Brune IYe = ye + s; 406d9c9ebe5SPeter Brune Ys = ys - s; 407d9c9ebe5SPeter Brune Ye = ye + s; 40888661749SPeter Brune } 40988661749SPeter Brune 410d9c9ebe5SPeter Brune if (zs-s > 0) { 411d9c9ebe5SPeter Brune Zs = zs - s; IZs = zs - s; 41288661749SPeter Brune } else { 4138865f1eaSKarl Rupp if (bz) Zs = zs - s; 4148865f1eaSKarl Rupp else Zs = 0; 41588661749SPeter Brune IZs = 0; 41688661749SPeter Brune } 417d9c9ebe5SPeter Brune if (ze+s <= P) { 418d9c9ebe5SPeter Brune Ze = ze + s; IZe = ze + s; 41988661749SPeter Brune } else { 4208865f1eaSKarl Rupp if (bz) Ze = ze + s; 4218865f1eaSKarl Rupp else Ze = P; 42288661749SPeter Brune IZe = P; 42388661749SPeter Brune } 42488661749SPeter Brune 425bff4a2f0SMatthew G. Knepley if (bz == DM_BOUNDARY_PERIODIC || bz == DM_BOUNDARY_MIRROR) { 426d9c9ebe5SPeter Brune IZs = zs - s; 427d9c9ebe5SPeter Brune IZe = ze + s; 428d9c9ebe5SPeter Brune Zs = zs - s; 429d9c9ebe5SPeter Brune Ze = ze + s; 43088661749SPeter Brune } 43147c6ae99SBarry Smith 43247c6ae99SBarry Smith /* Resize all X parameters to reflect w */ 433d9c9ebe5SPeter Brune s_x = s; 434d9c9ebe5SPeter Brune s_y = s; 435d9c9ebe5SPeter Brune s_z = s; 43647c6ae99SBarry Smith 43747c6ae99SBarry Smith /* determine starting point of each processor */ 43847c6ae99SBarry Smith nn = x*y*z; 4399566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(size+1,&bases,size,&ldims)); 4409566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm)); 44147c6ae99SBarry Smith bases[0] = 0; 4428865f1eaSKarl Rupp for (i=1; i<=size; i++) bases[i] = ldims[i-1]; 4438865f1eaSKarl Rupp for (i=1; i<=size; i++) bases[i] += bases[i-1]; 444ce00eea3SSatish Balay base = bases[rank]*dof; 44547c6ae99SBarry Smith 44647c6ae99SBarry Smith /* allocate the base parallel and sequential vectors */ 447ce00eea3SSatish Balay dd->Nlocal = x*y*z*dof; 4489566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,NULL,&global)); 449ce00eea3SSatish Balay dd->nlocal = (Xe-Xs)*(Ye-Ys)*(Ze-Zs)*dof; 4509566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,NULL,&local)); 45147c6ae99SBarry Smith 4524104a7a0SPatrick Sanan /* generate global to local vector scatter and local to global mapping*/ 45347c6ae99SBarry Smith 454ce00eea3SSatish Balay /* global to local must include ghost points within the domain, 455ce00eea3SSatish Balay but not ghost points outside the domain that aren't periodic */ 4569566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((IXe-IXs)*(IYe-IYs)*(IZe-IZs),&idx)); 457d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_BOX) { 458ce00eea3SSatish Balay left = IXs - Xs; right = left + (IXe-IXs); 459ce00eea3SSatish Balay bottom = IYs - Ys; top = bottom + (IYe-IYs); 460ce00eea3SSatish Balay down = IZs - Zs; up = down + (IZe-IZs); 461ce00eea3SSatish Balay count = 0; 462ce00eea3SSatish Balay for (i=down; i<up; i++) { 463ce00eea3SSatish Balay 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; 466ce00eea3SSatish Balay } 467ce00eea3SSatish Balay } 468ce00eea3SSatish Balay } 4699566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to)); 47047c6ae99SBarry Smith } else { 47147c6ae99SBarry Smith /* This is way ugly! We need to list the funny cross type region */ 472ce00eea3SSatish Balay left = xs - Xs; right = left + x; 47347c6ae99SBarry Smith bottom = ys - Ys; top = bottom + y; 47447c6ae99SBarry Smith down = zs - Zs; up = down + z; 47547c6ae99SBarry Smith count = 0; 476a5b23f4aSJose E. Roman /* the bottom chunk */ 477ce00eea3SSatish Balay for (i=(IZs-Zs); i<down; i++) { 47847c6ae99SBarry Smith for (j=bottom; j<top; j++) { 479ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 48047c6ae99SBarry Smith } 48147c6ae99SBarry Smith } 48247c6ae99SBarry Smith /* the middle piece */ 48347c6ae99SBarry Smith for (i=down; i<up; i++) { 48447c6ae99SBarry Smith /* front */ 485ce00eea3SSatish Balay for (j=(IYs-Ys); j<bottom; j++) { 486ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 48747c6ae99SBarry Smith } 48847c6ae99SBarry Smith /* middle */ 48947c6ae99SBarry Smith for (j=bottom; j<top; j++) { 490ce00eea3SSatish Balay for (k=IXs-Xs; k<IXe-Xs; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 49147c6ae99SBarry Smith } 49247c6ae99SBarry Smith /* back */ 493ce00eea3SSatish Balay for (j=top; j<top+IYe-ye; 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 top piece */ 498ce00eea3SSatish Balay for (i=up; i<up+IZe-ze; i++) { 49947c6ae99SBarry Smith for (j=bottom; j<top; j++) { 500ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 50147c6ae99SBarry Smith } 50247c6ae99SBarry Smith } 5039566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to)); 50447c6ae99SBarry Smith } 50547c6ae99SBarry Smith 50647c6ae99SBarry Smith /* determine who lies on each side of use stored in n24 n25 n26 50747c6ae99SBarry Smith n21 n22 n23 50847c6ae99SBarry Smith n18 n19 n20 50947c6ae99SBarry Smith 51047c6ae99SBarry Smith n15 n16 n17 51147c6ae99SBarry Smith n12 n14 51247c6ae99SBarry Smith n9 n10 n11 51347c6ae99SBarry Smith 51447c6ae99SBarry Smith n6 n7 n8 51547c6ae99SBarry Smith n3 n4 n5 51647c6ae99SBarry Smith n0 n1 n2 51747c6ae99SBarry Smith */ 51847c6ae99SBarry Smith 51947c6ae99SBarry Smith /* Solve for X,Y, and Z Periodic Case First, Then Modify Solution */ 52047c6ae99SBarry Smith /* Assume Nodes are Internal to the Cube */ 52147c6ae99SBarry Smith n0 = rank - m*n - m - 1; 52247c6ae99SBarry Smith n1 = rank - m*n - m; 52347c6ae99SBarry Smith n2 = rank - m*n - m + 1; 52447c6ae99SBarry Smith n3 = rank - m*n -1; 52547c6ae99SBarry Smith n4 = rank - m*n; 52647c6ae99SBarry Smith n5 = rank - m*n + 1; 52747c6ae99SBarry Smith n6 = rank - m*n + m - 1; 52847c6ae99SBarry Smith n7 = rank - m*n + m; 52947c6ae99SBarry Smith n8 = rank - m*n + m + 1; 53047c6ae99SBarry Smith 53147c6ae99SBarry Smith n9 = rank - m - 1; 53247c6ae99SBarry Smith n10 = rank - m; 53347c6ae99SBarry Smith n11 = rank - m + 1; 53447c6ae99SBarry Smith n12 = rank - 1; 53547c6ae99SBarry Smith n14 = rank + 1; 53647c6ae99SBarry Smith n15 = rank + m - 1; 53747c6ae99SBarry Smith n16 = rank + m; 53847c6ae99SBarry Smith n17 = rank + m + 1; 53947c6ae99SBarry Smith 54047c6ae99SBarry Smith n18 = rank + m*n - m - 1; 54147c6ae99SBarry Smith n19 = rank + m*n - m; 54247c6ae99SBarry Smith n20 = rank + m*n - m + 1; 54347c6ae99SBarry Smith n21 = rank + m*n - 1; 54447c6ae99SBarry Smith n22 = rank + m*n; 54547c6ae99SBarry Smith n23 = rank + m*n + 1; 54647c6ae99SBarry Smith n24 = rank + m*n + m - 1; 54747c6ae99SBarry Smith n25 = rank + m*n + m; 54847c6ae99SBarry Smith n26 = rank + m*n + m + 1; 54947c6ae99SBarry Smith 55047c6ae99SBarry Smith /* Assume Pieces are on Faces of Cube */ 55147c6ae99SBarry Smith 55247c6ae99SBarry Smith if (xs == 0) { /* First assume not corner or edge */ 55347c6ae99SBarry Smith n0 = rank -1 - (m*n); 55447c6ae99SBarry Smith n3 = rank + m -1 - (m*n); 55547c6ae99SBarry Smith n6 = rank + 2*m -1 - (m*n); 55647c6ae99SBarry Smith n9 = rank -1; 55747c6ae99SBarry Smith n12 = rank + m -1; 55847c6ae99SBarry Smith n15 = rank + 2*m -1; 55947c6ae99SBarry Smith n18 = rank -1 + (m*n); 56047c6ae99SBarry Smith n21 = rank + m -1 + (m*n); 56147c6ae99SBarry Smith n24 = rank + 2*m -1 + (m*n); 56247c6ae99SBarry Smith } 56347c6ae99SBarry Smith 564ce00eea3SSatish Balay if (xe == M) { /* First assume not corner or edge */ 56547c6ae99SBarry Smith n2 = rank -2*m +1 - (m*n); 56647c6ae99SBarry Smith n5 = rank - m +1 - (m*n); 56747c6ae99SBarry Smith n8 = rank +1 - (m*n); 56847c6ae99SBarry Smith n11 = rank -2*m +1; 56947c6ae99SBarry Smith n14 = rank - m +1; 57047c6ae99SBarry Smith n17 = rank +1; 57147c6ae99SBarry Smith n20 = rank -2*m +1 + (m*n); 57247c6ae99SBarry Smith n23 = rank - m +1 + (m*n); 57347c6ae99SBarry Smith n26 = rank +1 + (m*n); 57447c6ae99SBarry Smith } 57547c6ae99SBarry Smith 57647c6ae99SBarry Smith if (ys==0) { /* First assume not corner or edge */ 57747c6ae99SBarry Smith n0 = rank + m * (n-1) -1 - (m*n); 57847c6ae99SBarry Smith n1 = rank + m * (n-1) - (m*n); 57947c6ae99SBarry Smith n2 = rank + m * (n-1) +1 - (m*n); 58047c6ae99SBarry Smith n9 = rank + m * (n-1) -1; 58147c6ae99SBarry Smith n10 = rank + m * (n-1); 58247c6ae99SBarry Smith n11 = rank + m * (n-1) +1; 58347c6ae99SBarry Smith n18 = rank + m * (n-1) -1 + (m*n); 58447c6ae99SBarry Smith n19 = rank + m * (n-1) + (m*n); 58547c6ae99SBarry Smith n20 = rank + m * (n-1) +1 + (m*n); 58647c6ae99SBarry Smith } 58747c6ae99SBarry Smith 58847c6ae99SBarry Smith if (ye == N) { /* First assume not corner or edge */ 58947c6ae99SBarry Smith n6 = rank - m * (n-1) -1 - (m*n); 59047c6ae99SBarry Smith n7 = rank - m * (n-1) - (m*n); 59147c6ae99SBarry Smith n8 = rank - m * (n-1) +1 - (m*n); 59247c6ae99SBarry Smith n15 = rank - m * (n-1) -1; 59347c6ae99SBarry Smith n16 = rank - m * (n-1); 59447c6ae99SBarry Smith n17 = rank - m * (n-1) +1; 59547c6ae99SBarry Smith n24 = rank - m * (n-1) -1 + (m*n); 59647c6ae99SBarry Smith n25 = rank - m * (n-1) + (m*n); 59747c6ae99SBarry Smith n26 = rank - m * (n-1) +1 + (m*n); 59847c6ae99SBarry Smith } 59947c6ae99SBarry Smith 60047c6ae99SBarry Smith if (zs == 0) { /* First assume not corner or edge */ 60147c6ae99SBarry Smith n0 = size - (m*n) + rank - m - 1; 60247c6ae99SBarry Smith n1 = size - (m*n) + rank - m; 60347c6ae99SBarry Smith n2 = size - (m*n) + rank - m + 1; 60447c6ae99SBarry Smith n3 = size - (m*n) + rank - 1; 60547c6ae99SBarry Smith n4 = size - (m*n) + rank; 60647c6ae99SBarry Smith n5 = size - (m*n) + rank + 1; 60747c6ae99SBarry Smith n6 = size - (m*n) + rank + m - 1; 60847c6ae99SBarry Smith n7 = size - (m*n) + rank + m; 60947c6ae99SBarry Smith n8 = size - (m*n) + rank + m + 1; 61047c6ae99SBarry Smith } 61147c6ae99SBarry Smith 61247c6ae99SBarry Smith if (ze == P) { /* First assume not corner or edge */ 61347c6ae99SBarry Smith n18 = (m*n) - (size-rank) - m - 1; 61447c6ae99SBarry Smith n19 = (m*n) - (size-rank) - m; 61547c6ae99SBarry Smith n20 = (m*n) - (size-rank) - m + 1; 61647c6ae99SBarry Smith n21 = (m*n) - (size-rank) - 1; 61747c6ae99SBarry Smith n22 = (m*n) - (size-rank); 61847c6ae99SBarry Smith n23 = (m*n) - (size-rank) + 1; 61947c6ae99SBarry Smith n24 = (m*n) - (size-rank) + m - 1; 62047c6ae99SBarry Smith n25 = (m*n) - (size-rank) + m; 62147c6ae99SBarry Smith n26 = (m*n) - (size-rank) + m + 1; 62247c6ae99SBarry Smith } 62347c6ae99SBarry Smith 62447c6ae99SBarry Smith if ((xs==0) && (zs==0)) { /* Assume an edge, not corner */ 62547c6ae99SBarry Smith n0 = size - m*n + rank + m-1 - m; 62647c6ae99SBarry Smith n3 = size - m*n + rank + m-1; 62747c6ae99SBarry Smith n6 = size - m*n + rank + m-1 + m; 62847c6ae99SBarry Smith } 62947c6ae99SBarry Smith 63047c6ae99SBarry Smith if ((xs==0) && (ze==P)) { /* Assume an edge, not corner */ 63147c6ae99SBarry Smith n18 = m*n - (size - rank) + m-1 - m; 63247c6ae99SBarry Smith n21 = m*n - (size - rank) + m-1; 63347c6ae99SBarry Smith n24 = m*n - (size - rank) + m-1 + m; 63447c6ae99SBarry Smith } 63547c6ae99SBarry Smith 63647c6ae99SBarry Smith if ((xs==0) && (ys==0)) { /* Assume an edge, not corner */ 63747c6ae99SBarry Smith n0 = rank + m*n -1 - m*n; 63847c6ae99SBarry Smith n9 = rank + m*n -1; 63947c6ae99SBarry Smith n18 = rank + m*n -1 + m*n; 64047c6ae99SBarry Smith } 64147c6ae99SBarry Smith 64247c6ae99SBarry Smith if ((xs==0) && (ye==N)) { /* Assume an edge, not corner */ 64347c6ae99SBarry Smith n6 = rank - m*(n-1) + m-1 - m*n; 64447c6ae99SBarry Smith n15 = rank - m*(n-1) + m-1; 64547c6ae99SBarry Smith n24 = rank - m*(n-1) + m-1 + m*n; 64647c6ae99SBarry Smith } 64747c6ae99SBarry Smith 648ce00eea3SSatish Balay if ((xe==M) && (zs==0)) { /* Assume an edge, not corner */ 64947c6ae99SBarry Smith n2 = size - (m*n-rank) - (m-1) - m; 65047c6ae99SBarry Smith n5 = size - (m*n-rank) - (m-1); 65147c6ae99SBarry Smith n8 = size - (m*n-rank) - (m-1) + m; 65247c6ae99SBarry Smith } 65347c6ae99SBarry Smith 654ce00eea3SSatish Balay if ((xe==M) && (ze==P)) { /* Assume an edge, not corner */ 65547c6ae99SBarry Smith n20 = m*n - (size - rank) - (m-1) - m; 65647c6ae99SBarry Smith n23 = m*n - (size - rank) - (m-1); 65747c6ae99SBarry Smith n26 = m*n - (size - rank) - (m-1) + m; 65847c6ae99SBarry Smith } 65947c6ae99SBarry Smith 660ce00eea3SSatish Balay if ((xe==M) && (ys==0)) { /* Assume an edge, not corner */ 66147c6ae99SBarry Smith n2 = rank + m*(n-1) - (m-1) - m*n; 66247c6ae99SBarry Smith n11 = rank + m*(n-1) - (m-1); 66347c6ae99SBarry Smith n20 = rank + m*(n-1) - (m-1) + m*n; 66447c6ae99SBarry Smith } 66547c6ae99SBarry Smith 666ce00eea3SSatish Balay if ((xe==M) && (ye==N)) { /* Assume an edge, not corner */ 66747c6ae99SBarry Smith n8 = rank - m*n +1 - m*n; 66847c6ae99SBarry Smith n17 = rank - m*n +1; 66947c6ae99SBarry Smith n26 = rank - m*n +1 + m*n; 67047c6ae99SBarry Smith } 67147c6ae99SBarry Smith 67247c6ae99SBarry Smith if ((ys==0) && (zs==0)) { /* Assume an edge, not corner */ 67347c6ae99SBarry Smith n0 = size - m + rank -1; 67447c6ae99SBarry Smith n1 = size - m + rank; 67547c6ae99SBarry Smith n2 = size - m + rank +1; 67647c6ae99SBarry Smith } 67747c6ae99SBarry Smith 67847c6ae99SBarry Smith if ((ys==0) && (ze==P)) { /* Assume an edge, not corner */ 67947c6ae99SBarry Smith n18 = m*n - (size - rank) + m*(n-1) -1; 68047c6ae99SBarry Smith n19 = m*n - (size - rank) + m*(n-1); 68147c6ae99SBarry Smith n20 = m*n - (size - rank) + m*(n-1) +1; 68247c6ae99SBarry Smith } 68347c6ae99SBarry Smith 68447c6ae99SBarry Smith if ((ye==N) && (zs==0)) { /* Assume an edge, not corner */ 68547c6ae99SBarry Smith n6 = size - (m*n-rank) - m * (n-1) -1; 68647c6ae99SBarry Smith n7 = size - (m*n-rank) - m * (n-1); 68747c6ae99SBarry Smith n8 = size - (m*n-rank) - m * (n-1) +1; 68847c6ae99SBarry Smith } 68947c6ae99SBarry Smith 69047c6ae99SBarry Smith if ((ye==N) && (ze==P)) { /* Assume an edge, not corner */ 69147c6ae99SBarry Smith n24 = rank - (size-m) -1; 69247c6ae99SBarry Smith n25 = rank - (size-m); 69347c6ae99SBarry Smith n26 = rank - (size-m) +1; 69447c6ae99SBarry Smith } 69547c6ae99SBarry Smith 69647c6ae99SBarry Smith /* Check for Corners */ 6978865f1eaSKarl Rupp if ((xs==0) && (ys==0) && (zs==0)) n0 = size -1; 6988865f1eaSKarl Rupp if ((xs==0) && (ys==0) && (ze==P)) n18 = m*n-1; 6998865f1eaSKarl Rupp if ((xs==0) && (ye==N) && (zs==0)) n6 = (size-1)-m*(n-1); 7008865f1eaSKarl Rupp if ((xs==0) && (ye==N) && (ze==P)) n24 = m-1; 7018865f1eaSKarl Rupp if ((xe==M) && (ys==0) && (zs==0)) n2 = size-m; 7028865f1eaSKarl Rupp if ((xe==M) && (ys==0) && (ze==P)) n20 = m*n-m; 7038865f1eaSKarl Rupp if ((xe==M) && (ye==N) && (zs==0)) n8 = size-m*n; 7048865f1eaSKarl Rupp if ((xe==M) && (ye==N) && (ze==P)) n26 = 0; 70547c6ae99SBarry Smith 70647c6ae99SBarry Smith /* Check for when not X,Y, and Z Periodic */ 70747c6ae99SBarry Smith 70847c6ae99SBarry Smith /* If not X periodic */ 709bff4a2f0SMatthew G. Knepley if (bx != DM_BOUNDARY_PERIODIC) { 7108865f1eaSKarl Rupp if (xs==0) n0 = n3 = n6 = n9 = n12 = n15 = n18 = n21 = n24 = -2; 7118865f1eaSKarl Rupp if (xe==M) n2 = n5 = n8 = n11 = n14 = n17 = n20 = n23 = n26 = -2; 71247c6ae99SBarry Smith } 71347c6ae99SBarry Smith 71447c6ae99SBarry Smith /* If not Y periodic */ 715bff4a2f0SMatthew G. Knepley if (by != DM_BOUNDARY_PERIODIC) { 7168865f1eaSKarl Rupp if (ys==0) n0 = n1 = n2 = n9 = n10 = n11 = n18 = n19 = n20 = -2; 7178865f1eaSKarl Rupp if (ye==N) n6 = n7 = n8 = n15 = n16 = n17 = n24 = n25 = n26 = -2; 71847c6ae99SBarry Smith } 71947c6ae99SBarry Smith 72047c6ae99SBarry Smith /* If not Z periodic */ 721bff4a2f0SMatthew G. Knepley if (bz != DM_BOUNDARY_PERIODIC) { 7228865f1eaSKarl Rupp if (zs==0) n0 = n1 = n2 = n3 = n4 = n5 = n6 = n7 = n8 = -2; 7238865f1eaSKarl Rupp if (ze==P) n18 = n19 = n20 = n21 = n22 = n23 = n24 = n25 = n26 = -2; 72447c6ae99SBarry Smith } 72547c6ae99SBarry Smith 7269566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(27,&dd->neighbors)); 7278865f1eaSKarl Rupp 72847c6ae99SBarry Smith dd->neighbors[0] = n0; 72947c6ae99SBarry Smith dd->neighbors[1] = n1; 73047c6ae99SBarry Smith dd->neighbors[2] = n2; 73147c6ae99SBarry Smith dd->neighbors[3] = n3; 73247c6ae99SBarry Smith dd->neighbors[4] = n4; 73347c6ae99SBarry Smith dd->neighbors[5] = n5; 73447c6ae99SBarry Smith dd->neighbors[6] = n6; 73547c6ae99SBarry Smith dd->neighbors[7] = n7; 73647c6ae99SBarry Smith dd->neighbors[8] = n8; 73747c6ae99SBarry Smith dd->neighbors[9] = n9; 73847c6ae99SBarry Smith dd->neighbors[10] = n10; 73947c6ae99SBarry Smith dd->neighbors[11] = n11; 74047c6ae99SBarry Smith dd->neighbors[12] = n12; 74147c6ae99SBarry Smith dd->neighbors[13] = rank; 74247c6ae99SBarry Smith dd->neighbors[14] = n14; 74347c6ae99SBarry Smith dd->neighbors[15] = n15; 74447c6ae99SBarry Smith dd->neighbors[16] = n16; 74547c6ae99SBarry Smith dd->neighbors[17] = n17; 74647c6ae99SBarry Smith dd->neighbors[18] = n18; 74747c6ae99SBarry Smith dd->neighbors[19] = n19; 74847c6ae99SBarry Smith dd->neighbors[20] = n20; 74947c6ae99SBarry Smith dd->neighbors[21] = n21; 75047c6ae99SBarry Smith dd->neighbors[22] = n22; 75147c6ae99SBarry Smith dd->neighbors[23] = n23; 75247c6ae99SBarry Smith dd->neighbors[24] = n24; 75347c6ae99SBarry Smith dd->neighbors[25] = n25; 75447c6ae99SBarry Smith dd->neighbors[26] = n26; 75547c6ae99SBarry Smith 75647c6ae99SBarry Smith /* If star stencil then delete the corner neighbors */ 757d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_STAR) { 75847c6ae99SBarry Smith /* save information about corner neighbors */ 75947c6ae99SBarry Smith sn0 = n0; sn1 = n1; sn2 = n2; sn3 = n3; sn5 = n5; sn6 = n6; sn7 = n7; 76047c6ae99SBarry Smith sn8 = n8; sn9 = n9; sn11 = n11; sn15 = n15; sn17 = n17; sn18 = n18; 76147c6ae99SBarry Smith sn19 = n19; sn20 = n20; sn21 = n21; sn23 = n23; sn24 = n24; sn25 = n25; 76247c6ae99SBarry Smith sn26 = n26; 7638865f1eaSKarl Rupp n0 = n1 = n2 = n3 = n5 = n6 = n7 = n8 = n9 = n11 = n15 = n17 = n18 = n19 = n20 = n21 = n23 = n24 = n25 = n26 = -1; 76447c6ae99SBarry Smith } 76547c6ae99SBarry Smith 7669566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((Xe-Xs)*(Ye-Ys)*(Ze-Zs),&idx)); 76747c6ae99SBarry Smith 76847c6ae99SBarry Smith nn = 0; 76947c6ae99SBarry Smith /* Bottom Level */ 77047c6ae99SBarry Smith for (k=0; k<s_z; k++) { 77147c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 77247c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 773ce00eea3SSatish Balay x_t = lx[n0 % m]; 77447c6ae99SBarry Smith y_t = ly[(n0 % (m*n))/m]; 77547c6ae99SBarry Smith z_t = lz[n0 / (m*n)]; 77647c6ae99SBarry 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; 7778865f1eaSKarl 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 */ 7788865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 77947c6ae99SBarry Smith } 78047c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 78147c6ae99SBarry Smith x_t = x; 78247c6ae99SBarry Smith y_t = ly[(n1 % (m*n))/m]; 78347c6ae99SBarry Smith z_t = lz[n1 / (m*n)]; 78447c6ae99SBarry 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; 7858865f1eaSKarl Rupp if (twod && (s_t < 0)) s_t = bases[n1] + x_t*y_t*z_t - (s_y+1-i)*x_t; /* 2D case */ 7868865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 78747c6ae99SBarry Smith } 78847c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 789ce00eea3SSatish Balay x_t = lx[n2 % m]; 79047c6ae99SBarry Smith y_t = ly[(n2 % (m*n))/m]; 79147c6ae99SBarry Smith z_t = lz[n2 / (m*n)]; 79247c6ae99SBarry 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; 7938865f1eaSKarl Rupp if (twod && (s_t < 0)) s_t = bases[n2] + x_t*y_t*z_t - (s_y+1-i)*x_t; /* 2D case */ 7948865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 79547c6ae99SBarry Smith } 79647c6ae99SBarry Smith } 79747c6ae99SBarry Smith 79847c6ae99SBarry Smith for (i=0; i<y; i++) { 79947c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 800ce00eea3SSatish Balay x_t = lx[n3 % m]; 80147c6ae99SBarry Smith y_t = y; 80247c6ae99SBarry Smith z_t = lz[n3 / (m*n)]; 80347c6ae99SBarry Smith s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8048865f1eaSKarl 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 */ 8058865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 80647c6ae99SBarry Smith } 80747c6ae99SBarry Smith 80847c6ae99SBarry Smith if (n4 >= 0) { /* middle */ 80947c6ae99SBarry Smith x_t = x; 81047c6ae99SBarry Smith y_t = y; 81147c6ae99SBarry Smith z_t = lz[n4 / (m*n)]; 81247c6ae99SBarry Smith s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8138865f1eaSKarl 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 */ 8148865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 815bff4a2f0SMatthew G. Knepley } else if (bz == DM_BOUNDARY_MIRROR) { 816a6fd6b0aSBarry Smith for (j=0; j<x; j++) idx[nn++] = bases[rank] + j + x*i + (s_z - k - 1)*x*y; 81747c6ae99SBarry Smith } 81847c6ae99SBarry Smith 81947c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 820ce00eea3SSatish Balay x_t = lx[n5 % m]; 82147c6ae99SBarry Smith y_t = y; 82247c6ae99SBarry Smith z_t = lz[n5 / (m*n)]; 82347c6ae99SBarry Smith s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8248865f1eaSKarl 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 */ 8258865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 82647c6ae99SBarry Smith } 82747c6ae99SBarry Smith } 82847c6ae99SBarry Smith 82947c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 83047c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 831ce00eea3SSatish Balay x_t = lx[n6 % m]; 83247c6ae99SBarry Smith y_t = ly[(n6 % (m*n))/m]; 83347c6ae99SBarry Smith z_t = lz[n6 / (m*n)]; 83447c6ae99SBarry Smith s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8358865f1eaSKarl 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 */ 8368865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 83747c6ae99SBarry Smith } 83847c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 83947c6ae99SBarry Smith x_t = x; 84047c6ae99SBarry Smith y_t = ly[(n7 % (m*n))/m]; 84147c6ae99SBarry Smith z_t = lz[n7 / (m*n)]; 84247c6ae99SBarry Smith s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8438865f1eaSKarl 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 */ 8448865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 84547c6ae99SBarry Smith } 84647c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 847ce00eea3SSatish Balay x_t = lx[n8 % m]; 84847c6ae99SBarry Smith y_t = ly[(n8 % (m*n))/m]; 84947c6ae99SBarry Smith z_t = lz[n8 / (m*n)]; 85047c6ae99SBarry Smith s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8518865f1eaSKarl 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 */ 8528865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 85347c6ae99SBarry Smith } 85447c6ae99SBarry Smith } 85547c6ae99SBarry Smith } 85647c6ae99SBarry Smith 85747c6ae99SBarry Smith /* Middle Level */ 85847c6ae99SBarry Smith for (k=0; k<z; k++) { 85947c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 86047c6ae99SBarry Smith if (n9 >= 0) { /* left below */ 861ce00eea3SSatish Balay x_t = lx[n9 % m]; 86247c6ae99SBarry Smith y_t = ly[(n9 % (m*n))/m]; 86347c6ae99SBarry Smith /* z_t = z; */ 86447c6ae99SBarry Smith s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 8658865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 86647c6ae99SBarry Smith } 86747c6ae99SBarry Smith if (n10 >= 0) { /* directly below */ 86847c6ae99SBarry Smith x_t = x; 86947c6ae99SBarry Smith y_t = ly[(n10 % (m*n))/m]; 87047c6ae99SBarry Smith /* z_t = z; */ 87147c6ae99SBarry Smith s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 8728865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 873bff4a2f0SMatthew G. Knepley } else if (by == DM_BOUNDARY_MIRROR) { 874a6fd6b0aSBarry Smith for (j=0; j<x; j++) idx[nn++] = bases[rank] + j + k*x*y + (s_y - i)*x; 87547c6ae99SBarry Smith } 87647c6ae99SBarry Smith if (n11 >= 0) { /* right below */ 877ce00eea3SSatish Balay x_t = lx[n11 % m]; 87847c6ae99SBarry Smith y_t = ly[(n11 % (m*n))/m]; 87947c6ae99SBarry Smith /* z_t = z; */ 88047c6ae99SBarry Smith s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 8818865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 88247c6ae99SBarry Smith } 88347c6ae99SBarry Smith } 88447c6ae99SBarry Smith 88547c6ae99SBarry Smith for (i=0; i<y; i++) { 88647c6ae99SBarry Smith if (n12 >= 0) { /* directly left */ 887ce00eea3SSatish Balay x_t = lx[n12 % m]; 88847c6ae99SBarry Smith y_t = y; 88947c6ae99SBarry Smith /* z_t = z; */ 89047c6ae99SBarry Smith s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t; 8918865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 892bff4a2f0SMatthew G. Knepley } else if (bx == DM_BOUNDARY_MIRROR) { 893a6fd6b0aSBarry Smith for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + s_x - j - 1 + k*x*y + i*x; 89447c6ae99SBarry Smith } 89547c6ae99SBarry Smith 89647c6ae99SBarry Smith /* Interior */ 89747c6ae99SBarry Smith s_t = bases[rank] + i*x + k*x*y; 8988865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = s_t++; 89947c6ae99SBarry Smith 90047c6ae99SBarry Smith if (n14 >= 0) { /* directly right */ 901ce00eea3SSatish Balay x_t = lx[n14 % m]; 90247c6ae99SBarry Smith y_t = y; 90347c6ae99SBarry Smith /* z_t = z; */ 90447c6ae99SBarry Smith s_t = bases[n14] + i*x_t + k*x_t*y_t; 9058865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 906bff4a2f0SMatthew G. Knepley } else if (bx == DM_BOUNDARY_MIRROR) { 907a6fd6b0aSBarry Smith for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x - j - 1 + k*x*y + i*x; 90847c6ae99SBarry Smith } 90947c6ae99SBarry Smith } 91047c6ae99SBarry Smith 91147c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 91247c6ae99SBarry Smith if (n15 >= 0) { /* left above */ 913ce00eea3SSatish Balay x_t = lx[n15 % m]; 91447c6ae99SBarry Smith y_t = ly[(n15 % (m*n))/m]; 91547c6ae99SBarry Smith /* z_t = z; */ 91647c6ae99SBarry Smith s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t; 9178865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 91847c6ae99SBarry Smith } 91947c6ae99SBarry Smith if (n16 >= 0) { /* directly above */ 92047c6ae99SBarry Smith x_t = x; 92147c6ae99SBarry Smith y_t = ly[(n16 % (m*n))/m]; 92247c6ae99SBarry Smith /* z_t = z; */ 92347c6ae99SBarry Smith s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t; 9248865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 925bff4a2f0SMatthew G. Knepley } else if (by == DM_BOUNDARY_MIRROR) { 926a6fd6b0aSBarry Smith for (j=0; j<x; j++) idx[nn++] = bases[rank] + j + k*x*y + (y-i)*x; 92747c6ae99SBarry Smith } 92847c6ae99SBarry Smith if (n17 >= 0) { /* right above */ 929ce00eea3SSatish Balay x_t = lx[n17 % m]; 93047c6ae99SBarry Smith y_t = ly[(n17 % (m*n))/m]; 93147c6ae99SBarry Smith /* z_t = z; */ 93247c6ae99SBarry Smith s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t; 9338865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 93447c6ae99SBarry Smith } 93547c6ae99SBarry Smith } 93647c6ae99SBarry Smith } 93747c6ae99SBarry Smith 93847c6ae99SBarry Smith /* Upper Level */ 93947c6ae99SBarry Smith for (k=0; k<s_z; k++) { 94047c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 94147c6ae99SBarry Smith if (n18 >= 0) { /* left below */ 942ce00eea3SSatish Balay x_t = lx[n18 % m]; 94347c6ae99SBarry Smith y_t = ly[(n18 % (m*n))/m]; 94447c6ae99SBarry Smith /* z_t = lz[n18 / (m*n)]; */ 94547c6ae99SBarry Smith s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 9468865f1eaSKarl 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 */ 9478865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 94847c6ae99SBarry Smith } 94947c6ae99SBarry Smith if (n19 >= 0) { /* directly below */ 95047c6ae99SBarry Smith x_t = x; 95147c6ae99SBarry Smith y_t = ly[(n19 % (m*n))/m]; 95247c6ae99SBarry Smith /* z_t = lz[n19 / (m*n)]; */ 95347c6ae99SBarry Smith s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 9548865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n19] - (s_y+1-i)*x_t + x_t*y_t; /* 2d case */ 9558865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 95647c6ae99SBarry Smith } 95747c6ae99SBarry Smith if (n20 >= 0) { /* right below */ 958ce00eea3SSatish Balay x_t = lx[n20 % m]; 95947c6ae99SBarry Smith y_t = ly[(n20 % (m*n))/m]; 96047c6ae99SBarry Smith /* z_t = lz[n20 / (m*n)]; */ 96147c6ae99SBarry Smith s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 9628865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n20] - (s_y+1-i)*x_t + x_t*y_t; /* 2d case */ 9638865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 96447c6ae99SBarry Smith } 96547c6ae99SBarry Smith } 96647c6ae99SBarry Smith 96747c6ae99SBarry Smith for (i=0; i<y; i++) { 96847c6ae99SBarry Smith if (n21 >= 0) { /* directly left */ 969ce00eea3SSatish Balay x_t = lx[n21 % m]; 97047c6ae99SBarry Smith y_t = y; 97147c6ae99SBarry Smith /* z_t = lz[n21 / (m*n)]; */ 97247c6ae99SBarry Smith s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t; 9738865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n21] + (i+1)*x_t - s_x; /* 2d case */ 9748865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 97547c6ae99SBarry Smith } 97647c6ae99SBarry Smith 97747c6ae99SBarry Smith if (n22 >= 0) { /* middle */ 97847c6ae99SBarry Smith x_t = x; 97947c6ae99SBarry Smith y_t = y; 98047c6ae99SBarry Smith /* z_t = lz[n22 / (m*n)]; */ 98147c6ae99SBarry Smith s_t = bases[n22] + i*x_t + k*x_t*y_t; 9828865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n22] + i*x_t; /* 2d case */ 9838865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 984bff4a2f0SMatthew G. Knepley } else if (bz == DM_BOUNDARY_MIRROR) { 985a6fd6b0aSBarry Smith for (j=0; j<x; j++) idx[nn++] = bases[rank] + j + (z-k-1)*x*y + i*x; 98647c6ae99SBarry Smith } 98747c6ae99SBarry Smith 98847c6ae99SBarry Smith if (n23 >= 0) { /* directly right */ 989ce00eea3SSatish Balay x_t = lx[n23 % m]; 99047c6ae99SBarry Smith y_t = y; 99147c6ae99SBarry Smith /* z_t = lz[n23 / (m*n)]; */ 99247c6ae99SBarry Smith s_t = bases[n23] + i*x_t + k*x_t*y_t; 9938865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n23] + i*x_t; /* 2d case */ 9948865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 99547c6ae99SBarry Smith } 99647c6ae99SBarry Smith } 99747c6ae99SBarry Smith 99847c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 99947c6ae99SBarry Smith if (n24 >= 0) { /* left above */ 1000ce00eea3SSatish Balay x_t = lx[n24 % m]; 100147c6ae99SBarry Smith y_t = ly[(n24 % (m*n))/m]; 100247c6ae99SBarry Smith /* z_t = lz[n24 / (m*n)]; */ 100347c6ae99SBarry Smith s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t; 10048865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n24] + i*x_t - s_x; /* 2d case */ 10058865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 100647c6ae99SBarry Smith } 100747c6ae99SBarry Smith if (n25 >= 0) { /* directly above */ 100847c6ae99SBarry Smith x_t = x; 100947c6ae99SBarry Smith y_t = ly[(n25 % (m*n))/m]; 101047c6ae99SBarry Smith /* z_t = lz[n25 / (m*n)]; */ 101147c6ae99SBarry Smith s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t; 10128865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n25] + (i-1)*x_t; /* 2d case */ 10138865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 101447c6ae99SBarry Smith } 101547c6ae99SBarry Smith if (n26 >= 0) { /* right above */ 1016ce00eea3SSatish Balay x_t = lx[n26 % m]; 101747c6ae99SBarry Smith y_t = ly[(n26 % (m*n))/m]; 101847c6ae99SBarry Smith /* z_t = lz[n26 / (m*n)]; */ 101947c6ae99SBarry Smith s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t; 10208865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n26] + (i-1)*x_t; /* 2d case */ 10218865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 102247c6ae99SBarry Smith } 102347c6ae99SBarry Smith } 102447c6ae99SBarry Smith } 1025ce00eea3SSatish Balay 10269566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(comm,dof,nn,idx,PETSC_USE_POINTER,&from)); 10279566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(global,from,local,to,>ol)); 10289566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)da,(PetscObject)gtol)); 10299566063dSJacob Faibussowitsch PetscCall(ISDestroy(&to)); 10309566063dSJacob Faibussowitsch PetscCall(ISDestroy(&from)); 103147c6ae99SBarry Smith 1032d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_STAR) { 103347c6ae99SBarry Smith n0 = sn0; n1 = sn1; n2 = sn2; n3 = sn3; n5 = sn5; n6 = sn6; n7 = sn7; 103447c6ae99SBarry Smith n8 = sn8; n9 = sn9; n11 = sn11; n15 = sn15; n17 = sn17; n18 = sn18; 103547c6ae99SBarry Smith n19 = sn19; n20 = sn20; n21 = sn21; n23 = sn23; n24 = sn24; n25 = sn25; 103647c6ae99SBarry Smith n26 = sn26; 1037ce00eea3SSatish Balay } 103847c6ae99SBarry Smith 1039288e7d53SBarry Smith if (((stencil_type == DMDA_STENCIL_STAR) || (bx != DM_BOUNDARY_PERIODIC && bx) || (by != DM_BOUNDARY_PERIODIC && by) || (bz != DM_BOUNDARY_PERIODIC && bz))) { 1040ce00eea3SSatish Balay /* 1041ce00eea3SSatish Balay Recompute the local to global mappings, this time keeping the 1042ce00eea3SSatish Balay information about the cross corner processor numbers. 1043ce00eea3SSatish Balay */ 104447c6ae99SBarry Smith nn = 0; 104547c6ae99SBarry Smith /* Bottom Level */ 104647c6ae99SBarry Smith for (k=0; k<s_z; k++) { 104747c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 104847c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 1049ce00eea3SSatish Balay x_t = lx[n0 % m]; 105047c6ae99SBarry Smith y_t = ly[(n0 % (m*n))/m]; 105147c6ae99SBarry Smith z_t = lz[n0 / (m*n)]; 105247c6ae99SBarry 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; 10538865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1054ce00eea3SSatish Balay } else if (Xs-xs < 0 && Ys-ys < 0 && Zs-zs < 0) { 10558865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 105647c6ae99SBarry Smith } 105747c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 105847c6ae99SBarry Smith x_t = x; 105947c6ae99SBarry Smith y_t = ly[(n1 % (m*n))/m]; 106047c6ae99SBarry Smith z_t = lz[n1 / (m*n)]; 106147c6ae99SBarry 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; 10628865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1063ce00eea3SSatish Balay } else if (Ys-ys < 0 && Zs-zs < 0) { 10648865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 106547c6ae99SBarry Smith } 106647c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 1067ce00eea3SSatish Balay x_t = lx[n2 % m]; 106847c6ae99SBarry Smith y_t = ly[(n2 % (m*n))/m]; 106947c6ae99SBarry Smith z_t = lz[n2 / (m*n)]; 107047c6ae99SBarry 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; 10718865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1072ce00eea3SSatish Balay } else if (xe-Xe < 0 && Ys-ys < 0 && Zs-zs < 0) { 10738865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 107447c6ae99SBarry Smith } 107547c6ae99SBarry Smith } 107647c6ae99SBarry Smith 107747c6ae99SBarry Smith for (i=0; i<y; i++) { 107847c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 1079ce00eea3SSatish Balay x_t = lx[n3 % m]; 108047c6ae99SBarry Smith y_t = y; 108147c6ae99SBarry Smith z_t = lz[n3 / (m*n)]; 108247c6ae99SBarry Smith s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 10838865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1084ce00eea3SSatish Balay } else if (Xs-xs < 0 && Zs-zs < 0) { 10858865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 108647c6ae99SBarry Smith } 108747c6ae99SBarry Smith 108847c6ae99SBarry Smith if (n4 >= 0) { /* middle */ 108947c6ae99SBarry Smith x_t = x; 109047c6ae99SBarry Smith y_t = y; 109147c6ae99SBarry Smith z_t = lz[n4 / (m*n)]; 109247c6ae99SBarry Smith s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 10938865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1094ce00eea3SSatish Balay } else if (Zs-zs < 0) { 1095bff4a2f0SMatthew G. Knepley if (bz == DM_BOUNDARY_MIRROR) { 1096a6fd6b0aSBarry Smith for (j=0; j<x; j++) idx[nn++] = bases[rank] + j + x*i + (s_z - k - 1)*x*y; 1097c2859e5eSBarry Smith } else { 10988865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 109947c6ae99SBarry Smith } 1100c2859e5eSBarry Smith } 110147c6ae99SBarry Smith 110247c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 1103ce00eea3SSatish Balay x_t = lx[n5 % m]; 110447c6ae99SBarry Smith y_t = y; 110547c6ae99SBarry Smith z_t = lz[n5 / (m*n)]; 110647c6ae99SBarry Smith s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 11078865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1108ce00eea3SSatish Balay } else if (xe-Xe < 0 && Zs-zs < 0) { 11098865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 111047c6ae99SBarry Smith } 111147c6ae99SBarry Smith } 111247c6ae99SBarry Smith 111347c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 111447c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 1115ce00eea3SSatish Balay x_t = lx[n6 % m]; 111647c6ae99SBarry Smith y_t = ly[(n6 % (m*n))/m]; 111747c6ae99SBarry Smith z_t = lz[n6 / (m*n)]; 111847c6ae99SBarry Smith s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 11198865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1120ce00eea3SSatish Balay } else if (Xs-xs < 0 && ye-Ye < 0 && Zs-zs < 0) { 11218865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 112247c6ae99SBarry Smith } 112347c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 112447c6ae99SBarry Smith x_t = x; 112547c6ae99SBarry Smith y_t = ly[(n7 % (m*n))/m]; 112647c6ae99SBarry Smith z_t = lz[n7 / (m*n)]; 112747c6ae99SBarry Smith s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 11288865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1129ce00eea3SSatish Balay } else if (ye-Ye < 0 && Zs-zs < 0) { 11308865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 113147c6ae99SBarry Smith } 113247c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 1133ce00eea3SSatish Balay x_t = lx[n8 % m]; 113447c6ae99SBarry Smith y_t = ly[(n8 % (m*n))/m]; 113547c6ae99SBarry Smith z_t = lz[n8 / (m*n)]; 113647c6ae99SBarry Smith s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 11378865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1138ce00eea3SSatish Balay } else if (xe-Xe < 0 && ye-Ye < 0 && Zs-zs < 0) { 11398865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 114047c6ae99SBarry Smith } 114147c6ae99SBarry Smith } 114247c6ae99SBarry Smith } 114347c6ae99SBarry Smith 114447c6ae99SBarry Smith /* Middle Level */ 114547c6ae99SBarry Smith for (k=0; k<z; k++) { 114647c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 114747c6ae99SBarry Smith if (n9 >= 0) { /* left below */ 1148ce00eea3SSatish Balay x_t = lx[n9 % m]; 114947c6ae99SBarry Smith y_t = ly[(n9 % (m*n))/m]; 115047c6ae99SBarry Smith /* z_t = z; */ 115147c6ae99SBarry Smith s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 11528865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1153ce00eea3SSatish Balay } else if (Xs-xs < 0 && Ys-ys < 0) { 11548865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 115547c6ae99SBarry Smith } 115647c6ae99SBarry Smith if (n10 >= 0) { /* directly below */ 115747c6ae99SBarry Smith x_t = x; 115847c6ae99SBarry Smith y_t = ly[(n10 % (m*n))/m]; 115947c6ae99SBarry Smith /* z_t = z; */ 116047c6ae99SBarry Smith s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 11618865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1162ce00eea3SSatish Balay } else if (Ys-ys < 0) { 1163bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_MIRROR) { 1164feb237baSPierre Jolivet for (j=0; j<x; j++) idx[nn++] = bases[rank] + j + k*x*y + (s_y - i)*x; 1165c2859e5eSBarry Smith } else { 11668865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 1167c2859e5eSBarry Smith } 116847c6ae99SBarry Smith } 116947c6ae99SBarry Smith if (n11 >= 0) { /* right below */ 1170ce00eea3SSatish Balay x_t = lx[n11 % m]; 117147c6ae99SBarry Smith y_t = ly[(n11 % (m*n))/m]; 117247c6ae99SBarry Smith /* z_t = z; */ 117347c6ae99SBarry Smith s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 11748865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1175ce00eea3SSatish Balay } else if (xe-Xe < 0 && Ys-ys < 0) { 11768865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 117747c6ae99SBarry Smith } 117847c6ae99SBarry Smith } 117947c6ae99SBarry Smith 118047c6ae99SBarry Smith for (i=0; i<y; i++) { 118147c6ae99SBarry Smith if (n12 >= 0) { /* directly left */ 1182ce00eea3SSatish Balay x_t = lx[n12 % m]; 118347c6ae99SBarry Smith y_t = y; 118447c6ae99SBarry Smith /* z_t = z; */ 118547c6ae99SBarry Smith s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t; 11868865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1187ce00eea3SSatish Balay } else if (Xs-xs < 0) { 1188bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_MIRROR) { 1189a6fd6b0aSBarry Smith for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + s_x - j - 1 + k*x*y + i*x; 1190c2859e5eSBarry Smith } else { 11918865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 119247c6ae99SBarry Smith } 1193c2859e5eSBarry Smith } 119447c6ae99SBarry Smith 119547c6ae99SBarry Smith /* Interior */ 119647c6ae99SBarry Smith s_t = bases[rank] + i*x + k*x*y; 11978865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = s_t++; 119847c6ae99SBarry Smith 119947c6ae99SBarry Smith if (n14 >= 0) { /* directly right */ 1200ce00eea3SSatish Balay x_t = lx[n14 % m]; 120147c6ae99SBarry Smith y_t = y; 120247c6ae99SBarry Smith /* z_t = z; */ 120347c6ae99SBarry Smith s_t = bases[n14] + i*x_t + k*x_t*y_t; 12048865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1205ce00eea3SSatish Balay } else if (xe-Xe < 0) { 1206bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_MIRROR) { 1207a6fd6b0aSBarry Smith for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x - j - 1 + k*x*y + i*x; 1208c2859e5eSBarry Smith } else { 12098865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 121047c6ae99SBarry Smith } 121147c6ae99SBarry Smith } 1212c2859e5eSBarry Smith } 121347c6ae99SBarry Smith 121447c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 121547c6ae99SBarry Smith if (n15 >= 0) { /* left above */ 1216ce00eea3SSatish Balay x_t = lx[n15 % m]; 121747c6ae99SBarry Smith y_t = ly[(n15 % (m*n))/m]; 121847c6ae99SBarry Smith /* z_t = z; */ 121947c6ae99SBarry Smith s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t; 12208865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1221ce00eea3SSatish Balay } else if (Xs-xs < 0 && ye-Ye < 0) { 12228865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 122347c6ae99SBarry Smith } 122447c6ae99SBarry Smith if (n16 >= 0) { /* directly above */ 122547c6ae99SBarry Smith x_t = x; 122647c6ae99SBarry Smith y_t = ly[(n16 % (m*n))/m]; 122747c6ae99SBarry Smith /* z_t = z; */ 122847c6ae99SBarry Smith s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t; 12298865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1230ce00eea3SSatish Balay } else if (ye-Ye < 0) { 1231bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_MIRROR) { 1232a6fd6b0aSBarry Smith for (j=0; j<x; j++) idx[nn++] = bases[rank] + j + k*x*y + (y-i)*x; 1233c2859e5eSBarry Smith } else { 12348865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 123547c6ae99SBarry Smith } 1236c2859e5eSBarry Smith } 123747c6ae99SBarry Smith if (n17 >= 0) { /* right above */ 1238ce00eea3SSatish Balay x_t = lx[n17 % m]; 123947c6ae99SBarry Smith y_t = ly[(n17 % (m*n))/m]; 124047c6ae99SBarry Smith /* z_t = z; */ 124147c6ae99SBarry Smith s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t; 12428865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1243ce00eea3SSatish Balay } else if (xe-Xe < 0 && ye-Ye < 0) { 12448865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 124547c6ae99SBarry Smith } 124647c6ae99SBarry Smith } 124747c6ae99SBarry Smith } 124847c6ae99SBarry Smith 124947c6ae99SBarry Smith /* Upper Level */ 125047c6ae99SBarry Smith for (k=0; k<s_z; k++) { 125147c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 125247c6ae99SBarry Smith if (n18 >= 0) { /* left below */ 1253ce00eea3SSatish Balay x_t = lx[n18 % m]; 125447c6ae99SBarry Smith y_t = ly[(n18 % (m*n))/m]; 125547c6ae99SBarry Smith /* z_t = lz[n18 / (m*n)]; */ 125647c6ae99SBarry Smith s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 12578865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1258ce00eea3SSatish Balay } else if (Xs-xs < 0 && Ys-ys < 0 && ze-Ze < 0) { 12598865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 126047c6ae99SBarry Smith } 126147c6ae99SBarry Smith if (n19 >= 0) { /* directly below */ 126247c6ae99SBarry Smith x_t = x; 126347c6ae99SBarry Smith y_t = ly[(n19 % (m*n))/m]; 126447c6ae99SBarry Smith /* z_t = lz[n19 / (m*n)]; */ 126547c6ae99SBarry Smith s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 12668865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1267ce00eea3SSatish Balay } else if (Ys-ys < 0 && ze-Ze < 0) { 12688865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 126947c6ae99SBarry Smith } 127047c6ae99SBarry Smith if (n20 >= 0) { /* right below */ 1271ce00eea3SSatish Balay x_t = lx[n20 % m]; 127247c6ae99SBarry Smith y_t = ly[(n20 % (m*n))/m]; 127347c6ae99SBarry Smith /* z_t = lz[n20 / (m*n)]; */ 127447c6ae99SBarry Smith s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 12758865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1276ce00eea3SSatish Balay } else if (xe-Xe < 0 && Ys-ys < 0 && ze-Ze < 0) { 12778865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 127847c6ae99SBarry Smith } 127947c6ae99SBarry Smith } 128047c6ae99SBarry Smith 128147c6ae99SBarry Smith for (i=0; i<y; i++) { 128247c6ae99SBarry Smith if (n21 >= 0) { /* directly left */ 1283ce00eea3SSatish Balay x_t = lx[n21 % m]; 128447c6ae99SBarry Smith y_t = y; 128547c6ae99SBarry Smith /* z_t = lz[n21 / (m*n)]; */ 128647c6ae99SBarry Smith s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t; 12878865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1288ce00eea3SSatish Balay } else if (Xs-xs < 0 && ze-Ze < 0) { 12898865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 129047c6ae99SBarry Smith } 129147c6ae99SBarry Smith 129247c6ae99SBarry Smith if (n22 >= 0) { /* middle */ 129347c6ae99SBarry Smith x_t = x; 129447c6ae99SBarry Smith y_t = y; 129547c6ae99SBarry Smith /* z_t = lz[n22 / (m*n)]; */ 129647c6ae99SBarry Smith s_t = bases[n22] + i*x_t + k*x_t*y_t; 12978865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1298ce00eea3SSatish Balay } else if (ze-Ze < 0) { 1299bff4a2f0SMatthew G. Knepley if (bz == DM_BOUNDARY_MIRROR) { 1300a6fd6b0aSBarry Smith for (j=0; j<x; j++) idx[nn++] = bases[rank] + j + (z-k-1)*x*y + i*x; 1301c2859e5eSBarry Smith } else { 13028865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 130347c6ae99SBarry Smith } 1304c2859e5eSBarry Smith } 130547c6ae99SBarry Smith 130647c6ae99SBarry Smith if (n23 >= 0) { /* directly right */ 1307ce00eea3SSatish Balay x_t = lx[n23 % m]; 130847c6ae99SBarry Smith y_t = y; 130947c6ae99SBarry Smith /* z_t = lz[n23 / (m*n)]; */ 131047c6ae99SBarry Smith s_t = bases[n23] + i*x_t + k*x_t*y_t; 13118865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1312ce00eea3SSatish Balay } else if (xe-Xe < 0 && ze-Ze < 0) { 13138865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 131447c6ae99SBarry Smith } 131547c6ae99SBarry Smith } 131647c6ae99SBarry Smith 131747c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 131847c6ae99SBarry Smith if (n24 >= 0) { /* left above */ 1319ce00eea3SSatish Balay x_t = lx[n24 % m]; 132047c6ae99SBarry Smith y_t = ly[(n24 % (m*n))/m]; 132147c6ae99SBarry Smith /* z_t = lz[n24 / (m*n)]; */ 132247c6ae99SBarry Smith s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t; 13238865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1324ce00eea3SSatish Balay } else if (Xs-xs < 0 && ye-Ye < 0 && ze-Ze < 0) { 13258865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 132647c6ae99SBarry Smith } 132747c6ae99SBarry Smith if (n25 >= 0) { /* directly above */ 132847c6ae99SBarry Smith x_t = x; 132947c6ae99SBarry Smith y_t = ly[(n25 % (m*n))/m]; 133047c6ae99SBarry Smith /* z_t = lz[n25 / (m*n)]; */ 133147c6ae99SBarry Smith s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t; 13328865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1333ce00eea3SSatish Balay } else if (ye-Ye < 0 && ze-Ze < 0) { 13348865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 133547c6ae99SBarry Smith } 133647c6ae99SBarry Smith if (n26 >= 0) { /* right above */ 1337ce00eea3SSatish Balay x_t = lx[n26 % m]; 133847c6ae99SBarry Smith y_t = ly[(n26 % (m*n))/m]; 133947c6ae99SBarry Smith /* z_t = lz[n26 / (m*n)]; */ 134047c6ae99SBarry Smith s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t; 13418865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1342ce00eea3SSatish Balay } else if (xe-Xe < 0 && ye-Ye < 0 && ze-Ze < 0) { 13438865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 134447c6ae99SBarry Smith } 134547c6ae99SBarry Smith } 134647c6ae99SBarry Smith } 134747c6ae99SBarry Smith } 134847c6ae99SBarry Smith /* 134947c6ae99SBarry Smith Set the local to global ordering in the global vector, this allows use 135047c6ae99SBarry Smith of VecSetValuesLocal(). 135147c6ae99SBarry Smith */ 13529566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(comm,dof,nn,idx,PETSC_OWN_POINTER,&da->ltogmap)); 13539566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)da,(PetscObject)da->ltogmap)); 135447c6ae99SBarry Smith 13559566063dSJacob Faibussowitsch PetscCall(PetscFree2(bases,ldims)); 1356ce00eea3SSatish Balay dd->m = m; dd->n = n; dd->p = p; 1357ce00eea3SSatish Balay /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */ 1358ce00eea3SSatish Balay dd->xs = xs*dof; dd->xe = xe*dof; dd->ys = ys; dd->ye = ye; dd->zs = zs; dd->ze = ze; 1359ce00eea3SSatish Balay dd->Xs = Xs*dof; dd->Xe = Xe*dof; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = Zs; dd->Ze = Ze; 1360ce00eea3SSatish Balay 13619566063dSJacob Faibussowitsch PetscCall(VecDestroy(&local)); 13629566063dSJacob Faibussowitsch PetscCall(VecDestroy(&global)); 1363ce00eea3SSatish Balay 1364ce00eea3SSatish Balay dd->gtol = gtol; 1365ce00eea3SSatish Balay dd->base = base; 1366ce00eea3SSatish Balay da->ops->view = DMView_DA_3d; 13670298fd71SBarry Smith dd->ltol = NULL; 13680298fd71SBarry Smith dd->ao = NULL; 136947c6ae99SBarry Smith PetscFunctionReturn(0); 137047c6ae99SBarry Smith } 1371564755cdSBarry Smith 137247c6ae99SBarry Smith /*@C 1373aa219208SBarry Smith DMDACreate3d - Creates an object that will manage the communication of three-dimensional 137447c6ae99SBarry Smith regular array data that is distributed across some processors. 137547c6ae99SBarry Smith 1376d083f849SBarry Smith Collective 137747c6ae99SBarry Smith 137847c6ae99SBarry Smith Input Parameters: 137947c6ae99SBarry Smith + comm - MPI communicator 13801321219cSEthan Coon . bx,by,bz - type of ghost nodes the array have. 1381bff4a2f0SMatthew G. Knepley Use one of DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_PERIODIC. 1382aa219208SBarry Smith . stencil_type - Type of stencil (DMDA_STENCIL_STAR or DMDA_STENCIL_BOX) 1383897f7067SBarry Smith . M,N,P - global dimension in each direction of the array 138447c6ae99SBarry Smith . m,n,p - corresponding number of processors in each dimension 138547c6ae99SBarry Smith (or PETSC_DECIDE to have calculated) 138647c6ae99SBarry Smith . dof - number of degrees of freedom per node 138710d7c030SMatthew G Knepley . s - stencil width 138810d7c030SMatthew G Knepley - lx, ly, lz - arrays containing the number of nodes in each cell along 13890298fd71SBarry Smith the x, y, and z coordinates, or NULL. If non-null, these 139047c6ae99SBarry Smith must be of length as m,n,p and the corresponding 139147c6ae99SBarry Smith m,n, or p cannot be PETSC_DECIDE. Sum of the lx[] entries must be M, sum of 139247c6ae99SBarry Smith the ly[] must N, sum of the lz[] must be P 139347c6ae99SBarry Smith 139447c6ae99SBarry Smith Output Parameter: 139547c6ae99SBarry Smith . da - the resulting distributed array object 139647c6ae99SBarry Smith 139747c6ae99SBarry Smith Options Database Key: 1398706a11cbSBarry Smith + -dm_view - Calls DMView() at the conclusion of DMDACreate3d() 1399e43dc8daSBarry Smith . -da_grid_x <nx> - number of grid points in x direction 1400e43dc8daSBarry Smith . -da_grid_y <ny> - number of grid points in y direction 1401e43dc8daSBarry Smith . -da_grid_z <nz> - number of grid points in z direction 140247c6ae99SBarry Smith . -da_processors_x <MX> - number of processors in x direction 140347c6ae99SBarry Smith . -da_processors_y <MY> - number of processors in y direction 140447c6ae99SBarry Smith . -da_processors_z <MZ> - number of processors in z direction 1405e0f5d30fSBarry Smith . -da_refine_x <rx> - refinement ratio in x direction 1406e0f5d30fSBarry Smith . -da_refine_y <ry> - refinement ratio in y direction 1407e0f5d30fSBarry Smith . -da_refine_z <rz>- refinement ratio in z directio 1408e43dc8daSBarry Smith - -da_refine <n> - refine the DMDA n times before creating it 140947c6ae99SBarry Smith 141047c6ae99SBarry Smith Level: beginner 141147c6ae99SBarry Smith 141247c6ae99SBarry Smith Notes: 1413aa219208SBarry Smith The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the 1414aa219208SBarry Smith standard 7-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes 141547c6ae99SBarry Smith the standard 27-pt stencil. 141647c6ae99SBarry Smith 1417aa219208SBarry Smith The array data itself is NOT stored in the DMDA, it is stored in Vec objects; 1418564755cdSBarry Smith The appropriate vector objects can be obtained with calls to DMCreateGlobalVector() 1419564755cdSBarry Smith and DMCreateLocalVector() and calls to VecDuplicate() if more are needed. 142047c6ae99SBarry Smith 1421897f7067SBarry Smith You must call DMSetUp() after this call before using this DM. 1422897f7067SBarry Smith 1423897f7067SBarry Smith If you wish to use the options database to change values in the DMDA call DMSetFromOptions() after this call 1424897f7067SBarry Smith but before DMSetUp(). 1425897f7067SBarry Smith 1426aa219208SBarry Smith .seealso: DMDestroy(), DMView(), DMDACreate1d(), DMDACreate2d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(), 142799f0b487SRichard Tran Mills DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMLocalToLocalBegin(), DMLocalToLocalEnd(), DMDASetRefinementFactor(), 14283bd220d7SPatrick Sanan DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges(), 14293bd220d7SPatrick Sanan DMStagCreate3d() 143047c6ae99SBarry Smith 143147c6ae99SBarry Smith @*/ 1432bff4a2f0SMatthew G. Knepley PetscErrorCode DMDACreate3d(MPI_Comm comm,DMBoundaryType bx,DMBoundaryType by,DMBoundaryType bz,DMDAStencilType stencil_type,PetscInt M, 14339a42bb27SBarry 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) 143447c6ae99SBarry Smith { 143547c6ae99SBarry Smith PetscFunctionBegin; 14369566063dSJacob Faibussowitsch PetscCall(DMDACreate(comm, da)); 14379566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*da, 3)); 14389566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(*da, M, N, P)); 14399566063dSJacob Faibussowitsch PetscCall(DMDASetNumProcs(*da, m, n, p)); 14409566063dSJacob Faibussowitsch PetscCall(DMDASetBoundaryType(*da, bx, by, bz)); 14419566063dSJacob Faibussowitsch PetscCall(DMDASetDof(*da, dof)); 14429566063dSJacob Faibussowitsch PetscCall(DMDASetStencilType(*da, stencil_type)); 14439566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(*da, s)); 14449566063dSJacob Faibussowitsch PetscCall(DMDASetOwnershipRanges(*da, lx, ly, lz)); 144547c6ae99SBarry Smith PetscFunctionReturn(0); 144647c6ae99SBarry Smith } 1447