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; 5063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Load Balance - Grid Points: Min %" PetscInt_FMT " avg %" PetscInt_FMT " max %" PetscInt_FMT "\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)); 5663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"Processor [%d] M %" PetscInt_FMT " N %" PetscInt_FMT " P %" PetscInt_FMT " m %" PetscInt_FMT " n %" PetscInt_FMT " p %" PetscInt_FMT " w %" PetscInt_FMT " s %" PetscInt_FMT "\n",rank,dd->M,dd->N,dd->P,dd->m,dd->n,dd->p,dd->w,dd->s)); 5763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"X range of indices: %" PetscInt_FMT " %" PetscInt_FMT ", Y range of indices: %" PetscInt_FMT " %" PetscInt_FMT ", Z range of indices: %" PetscInt_FMT " %" PetscInt_FMT "\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)); 72*1baa6e33SBarry Smith } else if (format == PETSC_VIEWER_ASCII_GLVIS) PetscCall(DMView_DA_GLVis(da,viewer)); 73*1baa6e33SBarry Smith else PetscCall(DMView_DA_VTK(da,viewer)); 7447c6ae99SBarry Smith } else if (isdraw) { 7547c6ae99SBarry Smith PetscDraw draw; 7647c6ae99SBarry Smith PetscReal ymin = -1.0,ymax = (PetscReal)dd->N; 7747c6ae99SBarry Smith PetscReal xmin = -1.0,xmax = (PetscReal)((dd->M+2)*dd->P),x,y,ycoord,xcoord; 788ea3bf28SBarry Smith PetscInt k,plane,base; 798ea3bf28SBarry Smith const PetscInt *idx; 8047c6ae99SBarry Smith char node[10]; 8147c6ae99SBarry Smith PetscBool isnull; 8247c6ae99SBarry Smith 839566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(viewer,0,&draw)); 849566063dSJacob Faibussowitsch PetscCall(PetscDrawIsNull(draw,&isnull)); 855b399a63SLisandro Dalcin if (isnull) PetscFunctionReturn(0); 8647c6ae99SBarry Smith 879566063dSJacob Faibussowitsch PetscCall(PetscDrawCheckResizedWindow(draw)); 889566063dSJacob Faibussowitsch PetscCall(PetscDrawClear(draw)); 899566063dSJacob Faibussowitsch PetscCall(PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax)); 905b399a63SLisandro Dalcin 91d0609cedSBarry Smith PetscDrawCollectiveBegin(draw); 9247c6ae99SBarry Smith /* first processor draw all node lines */ 93dd400576SPatrick Sanan if (rank == 0) { 9447c6ae99SBarry Smith for (k=0; k<dd->P; k++) { 9547c6ae99SBarry Smith ymin = 0.0; ymax = (PetscReal)(dd->N - 1); 9647c6ae99SBarry Smith for (xmin=(PetscReal)(k*(dd->M+1)); xmin<(PetscReal)(dd->M+(k*(dd->M+1))); xmin++) { 979566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK)); 9847c6ae99SBarry Smith } 9947c6ae99SBarry Smith xmin = (PetscReal)(k*(dd->M+1)); xmax = xmin + (PetscReal)(dd->M - 1); 10047c6ae99SBarry Smith for (ymin=0; ymin<(PetscReal)dd->N; ymin++) { 1019566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK)); 10247c6ae99SBarry Smith } 10347c6ae99SBarry Smith } 10447c6ae99SBarry Smith } 105d0609cedSBarry Smith PetscDrawCollectiveEnd(draw); 1069566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 1079566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw)); 10847c6ae99SBarry Smith 109d0609cedSBarry Smith PetscDrawCollectiveBegin(draw); 1105b399a63SLisandro Dalcin /*Go through and draw for each plane*/ 1115b399a63SLisandro Dalcin for (k=0; k<dd->P; k++) { 11247c6ae99SBarry Smith if ((k >= dd->zs) && (k < dd->ze)) { 11347c6ae99SBarry Smith /* draw my box */ 11447c6ae99SBarry Smith ymin = dd->ys; 11547c6ae99SBarry Smith ymax = dd->ye - 1; 11647c6ae99SBarry Smith xmin = dd->xs/dd->w + (dd->M+1)*k; 11747c6ae99SBarry Smith xmax =(dd->xe-1)/dd->w + (dd->M+1)*k; 11847c6ae99SBarry Smith 1199566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED)); 1209566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED)); 1219566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED)); 1229566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED)); 12347c6ae99SBarry Smith 12447c6ae99SBarry Smith xmin = dd->xs/dd->w; 12547c6ae99SBarry Smith xmax =(dd->xe-1)/dd->w; 12647c6ae99SBarry Smith 127832b7cebSLisandro Dalcin /* identify which processor owns the box */ 1289566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(node,sizeof(node),"%d",(int)rank)); 1299566063dSJacob Faibussowitsch PetscCall(PetscDrawString(draw,xmin+(dd->M+1)*k+.2,ymin+.3,PETSC_DRAW_RED,node)); 13047c6ae99SBarry Smith /* put in numbers*/ 13147c6ae99SBarry Smith base = (dd->base+(dd->xe-dd->xs)*(dd->ye-dd->ys)*(k-dd->zs))/dd->w; 13247c6ae99SBarry Smith for (y=ymin; y<=ymax; y++) { 13347c6ae99SBarry Smith for (x=xmin+(dd->M+1)*k; x<=xmax+(dd->M+1)*k; x++) { 1349566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(node,sizeof(node),"%d",(int)base++)); 1359566063dSJacob Faibussowitsch PetscCall(PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node)); 13647c6ae99SBarry Smith } 13747c6ae99SBarry Smith } 13847c6ae99SBarry Smith 13947c6ae99SBarry Smith } 14047c6ae99SBarry Smith } 141d0609cedSBarry Smith PetscDrawCollectiveEnd(draw); 1429566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 1439566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw)); 14447c6ae99SBarry Smith 145d0609cedSBarry Smith PetscDrawCollectiveBegin(draw); 14647c6ae99SBarry Smith for (k=0-dd->s; k<dd->P+dd->s; k++) { 14747c6ae99SBarry Smith /* Go through and draw for each plane */ 14847c6ae99SBarry Smith if ((k >= dd->Zs) && (k < dd->Ze)) { 14947c6ae99SBarry Smith /* overlay ghost numbers, useful for error checking */ 150565245c5SBarry Smith base = (dd->Xe-dd->Xs)*(dd->Ye-dd->Ys)*(k-dd->Zs)/dd->w; 1519566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(da->ltogmap,&idx)); 15247c6ae99SBarry Smith plane=k; 153565245c5SBarry Smith /* Keep z wrap around points on the drawing */ 1548865f1eaSKarl Rupp if (k<0) plane=dd->P+k; 1558865f1eaSKarl Rupp if (k>=dd->P) plane=k-dd->P; 15647c6ae99SBarry Smith ymin = dd->Ys; ymax = dd->Ye; 15747c6ae99SBarry Smith xmin = (dd->M+1)*plane*dd->w; 15847c6ae99SBarry Smith xmax = (dd->M+1)*plane*dd->w+dd->M*dd->w; 15947c6ae99SBarry Smith for (y=ymin; y<ymax; y++) { 16047c6ae99SBarry Smith for (x=xmin+dd->Xs; x<xmin+dd->Xe; x+=dd->w) { 1618ea3bf28SBarry Smith sprintf(node,"%d",(int)(idx[base])); 16247c6ae99SBarry Smith ycoord = y; 16347c6ae99SBarry Smith /*Keep y wrap around points on drawing */ 1648865f1eaSKarl Rupp if (y<0) ycoord = dd->N+y; 1658865f1eaSKarl Rupp if (y>=dd->N) ycoord = y-dd->N; 16647c6ae99SBarry Smith xcoord = x; /* Keep x wrap points on drawing */ 1678865f1eaSKarl Rupp if (x<xmin) xcoord = xmax - (xmin-x); 1688865f1eaSKarl Rupp if (x>=xmax) xcoord = xmin + (x-xmax); 1699566063dSJacob Faibussowitsch PetscCall(PetscDrawString(draw,xcoord/dd->w,ycoord,PETSC_DRAW_BLUE,node)); 170565245c5SBarry Smith base++; 17147c6ae99SBarry Smith } 17247c6ae99SBarry Smith } 1739566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(da->ltogmap,&idx)); 17447c6ae99SBarry Smith } 17547c6ae99SBarry Smith } 176d0609cedSBarry Smith PetscDrawCollectiveEnd(draw); 1779566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 1789566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw)); 1799566063dSJacob Faibussowitsch PetscCall(PetscDrawSave(draw)); 180c9493c35SStefano Zampini } else if (isglvis) { 1819566063dSJacob Faibussowitsch PetscCall(DMView_DA_GLVis(da,viewer)); 1829a42bb27SBarry Smith } else if (isbinary) { 1839566063dSJacob Faibussowitsch PetscCall(DMView_DA_Binary(da,viewer)); 1849a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 1859a42bb27SBarry Smith } else if (ismatlab) { 1869566063dSJacob Faibussowitsch PetscCall(DMView_DA_Matlab(da,viewer)); 1879a42bb27SBarry Smith #endif 18811aeaf0aSBarry Smith } 18947c6ae99SBarry Smith PetscFunctionReturn(0); 19047c6ae99SBarry Smith } 19147c6ae99SBarry Smith 1927087cfbeSBarry Smith PetscErrorCode DMSetUp_DA_3D(DM da) 19347c6ae99SBarry Smith { 19447c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 19547c6ae99SBarry Smith const PetscInt M = dd->M; 19647c6ae99SBarry Smith const PetscInt N = dd->N; 19747c6ae99SBarry Smith const PetscInt P = dd->P; 19847c6ae99SBarry Smith PetscInt m = dd->m; 19947c6ae99SBarry Smith PetscInt n = dd->n; 20047c6ae99SBarry Smith PetscInt p = dd->p; 20147c6ae99SBarry Smith const PetscInt dof = dd->w; 20247c6ae99SBarry Smith const PetscInt s = dd->s; 203bff4a2f0SMatthew G. Knepley DMBoundaryType bx = dd->bx; 204bff4a2f0SMatthew G. Knepley DMBoundaryType by = dd->by; 205bff4a2f0SMatthew G. Knepley DMBoundaryType bz = dd->bz; 20619fd82e9SBarry Smith DMDAStencilType stencil_type = dd->stencil_type; 20747c6ae99SBarry Smith PetscInt *lx = dd->lx; 20847c6ae99SBarry Smith PetscInt *ly = dd->ly; 20947c6ae99SBarry Smith PetscInt *lz = dd->lz; 21047c6ae99SBarry Smith MPI_Comm comm; 21147c6ae99SBarry Smith PetscMPIInt rank,size; 212ce00eea3SSatish Balay PetscInt xs = 0,xe,ys = 0,ye,zs = 0,ze,x = 0,y = 0,z = 0; 213bd1fc5aeSBarry Smith PetscInt Xs,Xe,Ys,Ye,Zs,Ze,IXs,IXe,IYs,IYe,IZs,IZe,pm; 2148ea3bf28SBarry Smith PetscInt left,right,up,down,bottom,top,i,j,k,*idx,nn; 21547c6ae99SBarry Smith PetscInt n0,n1,n2,n3,n4,n5,n6,n7,n8,n9,n10,n11,n12,n14; 21647c6ae99SBarry Smith PetscInt n15,n16,n17,n18,n19,n20,n21,n22,n23,n24,n25,n26; 217db87c5ecSEthan Coon PetscInt *bases,*ldims,base,x_t,y_t,z_t,s_t,count,s_x,s_y,s_z; 21847c6ae99SBarry Smith PetscInt sn0 = 0,sn1 = 0,sn2 = 0,sn3 = 0,sn5 = 0,sn6 = 0,sn7 = 0; 21947c6ae99SBarry Smith PetscInt sn8 = 0,sn9 = 0,sn11 = 0,sn15 = 0,sn24 = 0,sn25 = 0,sn26 = 0; 22047c6ae99SBarry Smith PetscInt sn17 = 0,sn18 = 0,sn19 = 0,sn20 = 0,sn21 = 0,sn23 = 0; 22147c6ae99SBarry Smith Vec local,global; 222bd1fc5aeSBarry Smith VecScatter gtol; 22345b6f7e9SBarry Smith IS to,from; 2246f951b95Secoon PetscBool twod; 22547c6ae99SBarry Smith 22647c6ae99SBarry Smith PetscFunctionBegin; 2277a8be351SBarry 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"); 2289566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject) da, &comm)); 2293855c12bSBarry Smith #if !defined(PETSC_USE_64BIT_INDICES) 23063a3b9bcSJacob Faibussowitsch PetscCheck(((PetscInt64) M)*((PetscInt64) N)*((PetscInt64) P)*((PetscInt64) dof) <= (PetscInt64) PETSC_MPI_INT_MAX,comm,PETSC_ERR_INT_OVERFLOW,"Mesh of %" PetscInt_FMT " by %" PetscInt_FMT " by %" PetscInt_FMT " by %" PetscInt_FMT " (dof) is too large for 32 bit indices",M,N,P,dof); 2313855c12bSBarry Smith #endif 2323855c12bSBarry Smith 2339566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 2349566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm,&rank)); 23547c6ae99SBarry Smith 23647c6ae99SBarry Smith if (m != PETSC_DECIDE) { 23763a3b9bcSJacob Faibussowitsch PetscCheck(m >= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %" PetscInt_FMT,m); 23863a3b9bcSJacob Faibussowitsch else PetscCheck(m <= size,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %" PetscInt_FMT " %d",m,size); 23947c6ae99SBarry Smith } 24047c6ae99SBarry Smith if (n != PETSC_DECIDE) { 24163a3b9bcSJacob Faibussowitsch PetscCheck(n >= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %" PetscInt_FMT,n); 24263a3b9bcSJacob Faibussowitsch else PetscCheck(n <= size,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %" PetscInt_FMT " %d",n,size); 24347c6ae99SBarry Smith } 24447c6ae99SBarry Smith if (p != PETSC_DECIDE) { 24563a3b9bcSJacob Faibussowitsch PetscCheck(p >= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Z direction: %" PetscInt_FMT,p); 24663a3b9bcSJacob Faibussowitsch else PetscCheck(p <= size,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Z direction: %" PetscInt_FMT " %d",p,size); 24747c6ae99SBarry Smith } 2481dca8a05SBarry Smith PetscCheck(m <= 0 || n <= 0 || p <= 0 || m*n*p == size,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"m %" PetscInt_FMT " * n %" PetscInt_FMT " * p %" PetscInt_FMT " != size %d",m,n,p,size); 24947c6ae99SBarry Smith 25047c6ae99SBarry Smith /* Partition the array among the processors */ 25147c6ae99SBarry Smith if (m == PETSC_DECIDE && n != PETSC_DECIDE && p != PETSC_DECIDE) { 25247c6ae99SBarry Smith m = size/(n*p); 25347c6ae99SBarry Smith } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) { 25447c6ae99SBarry Smith n = size/(m*p); 25547c6ae99SBarry Smith } else if (m != PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) { 25647c6ae99SBarry Smith p = size/(m*n); 25747c6ae99SBarry Smith } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) { 25847c6ae99SBarry Smith /* try for squarish distribution */ 259369cc0aeSBarry Smith m = (int)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)N*p))); 26047c6ae99SBarry Smith if (!m) m = 1; 26147c6ae99SBarry Smith while (m > 0) { 26247c6ae99SBarry Smith n = size/(m*p); 26347c6ae99SBarry Smith if (m*n*p == size) break; 26447c6ae99SBarry Smith m--; 26547c6ae99SBarry Smith } 26663a3b9bcSJacob Faibussowitsch PetscCheck(m,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad p value: p = %" PetscInt_FMT,p); 26747c6ae99SBarry Smith if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;} 26847c6ae99SBarry Smith } else if (m == PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) { 26947c6ae99SBarry Smith /* try for squarish distribution */ 270369cc0aeSBarry Smith m = (int)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)P*n))); 27147c6ae99SBarry Smith if (!m) m = 1; 27247c6ae99SBarry Smith while (m > 0) { 27347c6ae99SBarry Smith p = size/(m*n); 27447c6ae99SBarry Smith if (m*n*p == size) break; 27547c6ae99SBarry Smith m--; 27647c6ae99SBarry Smith } 27763a3b9bcSJacob Faibussowitsch PetscCheck(m,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad n value: n = %" PetscInt_FMT,n); 27847c6ae99SBarry Smith if (M > P && m < p) {PetscInt _m = m; m = p; p = _m;} 27947c6ae99SBarry Smith } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) { 28047c6ae99SBarry Smith /* try for squarish distribution */ 281369cc0aeSBarry Smith n = (int)(0.5 + PetscSqrtReal(((PetscReal)N)*((PetscReal)size)/((PetscReal)P*m))); 28247c6ae99SBarry Smith if (!n) n = 1; 28347c6ae99SBarry Smith while (n > 0) { 28447c6ae99SBarry Smith p = size/(m*n); 28547c6ae99SBarry Smith if (m*n*p == size) break; 28647c6ae99SBarry Smith n--; 28747c6ae99SBarry Smith } 28863a3b9bcSJacob Faibussowitsch PetscCheck(n,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad m value: m = %" PetscInt_FMT,n); 28947c6ae99SBarry Smith if (N > P && n < p) {PetscInt _n = n; n = p; p = _n;} 29047c6ae99SBarry Smith } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) { 29147c6ae99SBarry Smith /* try for squarish distribution */ 292369cc0aeSBarry Smith n = (PetscInt)(0.5 + PetscPowReal(((PetscReal)N*N)*((PetscReal)size)/((PetscReal)P*M),(PetscReal)(1./3.))); 29347c6ae99SBarry Smith if (!n) n = 1; 29447c6ae99SBarry Smith while (n > 0) { 29547c6ae99SBarry Smith pm = size/n; 29647c6ae99SBarry Smith if (n*pm == size) break; 29747c6ae99SBarry Smith n--; 29847c6ae99SBarry Smith } 29947c6ae99SBarry Smith if (!n) n = 1; 300369cc0aeSBarry Smith m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)P*n))); 30147c6ae99SBarry Smith if (!m) m = 1; 30247c6ae99SBarry Smith while (m > 0) { 30347c6ae99SBarry Smith p = size/(m*n); 30447c6ae99SBarry Smith if (m*n*p == size) break; 30547c6ae99SBarry Smith m--; 30647c6ae99SBarry Smith } 30747c6ae99SBarry Smith if (M > P && m < p) {PetscInt _m = m; m = p; p = _m;} 30808401ef6SPierre Jolivet } else PetscCheck(m*n*p == size,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition"); 30947c6ae99SBarry Smith 31008401ef6SPierre Jolivet PetscCheck(m*n*p == size,PetscObjectComm((PetscObject)da),PETSC_ERR_PLIB,"Could not find good partition"); 31163a3b9bcSJacob Faibussowitsch PetscCheck(M >= m,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %" PetscInt_FMT " %" PetscInt_FMT,M,m); 31263a3b9bcSJacob Faibussowitsch PetscCheck(N >= n,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %" PetscInt_FMT " %" PetscInt_FMT,N,n); 31363a3b9bcSJacob Faibussowitsch PetscCheck(P >= p,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Partition in z direction is too fine! %" PetscInt_FMT " %" PetscInt_FMT,P,p); 31447c6ae99SBarry Smith 31547c6ae99SBarry Smith /* 31647c6ae99SBarry Smith Determine locally owned region 31747c6ae99SBarry Smith [x, y, or z]s is the first local node number, [x, y, z] is the number of local nodes 31847c6ae99SBarry Smith */ 31947c6ae99SBarry Smith 32047c6ae99SBarry Smith if (!lx) { 3219566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m, &dd->lx)); 32247c6ae99SBarry Smith lx = dd->lx; 3238865f1eaSKarl Rupp for (i=0; i<m; i++) lx[i] = M/m + ((M % m) > (i % m)); 32447c6ae99SBarry Smith } 32547c6ae99SBarry Smith x = lx[rank % m]; 32647c6ae99SBarry Smith xs = 0; 3278865f1eaSKarl Rupp for (i=0; i<(rank%m); i++) xs += lx[i]; 3281dca8a05SBarry Smith PetscCheck(x >= s || (m <= 1 && bx != DM_BOUNDARY_PERIODIC),PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local x-width of domain x %" PetscInt_FMT " is smaller than stencil width s %" PetscInt_FMT,x,s); 32947c6ae99SBarry Smith 33047c6ae99SBarry Smith if (!ly) { 3319566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &dd->ly)); 33247c6ae99SBarry Smith ly = dd->ly; 3338865f1eaSKarl Rupp for (i=0; i<n; i++) ly[i] = N/n + ((N % n) > (i % n)); 33447c6ae99SBarry Smith } 33547c6ae99SBarry Smith y = ly[(rank % (m*n))/m]; 3361dca8a05SBarry Smith PetscCheck(y >= s || (n <= 1 && by != DM_BOUNDARY_PERIODIC),PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local y-width of domain y %" PetscInt_FMT " is smaller than stencil width s %" PetscInt_FMT,y,s); 33730729d88SBarry Smith 33847c6ae99SBarry Smith ys = 0; 3398865f1eaSKarl Rupp for (i=0; i<(rank % (m*n))/m; i++) ys += ly[i]; 34047c6ae99SBarry Smith 34147c6ae99SBarry Smith if (!lz) { 3429566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(p, &dd->lz)); 34347c6ae99SBarry Smith lz = dd->lz; 3448865f1eaSKarl Rupp for (i=0; i<p; i++) lz[i] = P/p + ((P % p) > (i % p)); 34547c6ae99SBarry Smith } 34647c6ae99SBarry Smith z = lz[rank/(m*n)]; 347bcea557cSEthan Coon 348fdc81ce8SEthan Coon /* note this is different than x- and y-, as we will handle as an important special 349bff4a2f0SMatthew G. Knepley case when p=P=1 and DM_BOUNDARY_PERIODIC and s > z. This is to deal with 2D problems 350fdc81ce8SEthan Coon in a 3D code. Additional code for this case is noted with "2d case" comments */ 3516f951b95Secoon twod = PETSC_FALSE; 3528865f1eaSKarl Rupp if (P == 1) twod = PETSC_TRUE; 3531dca8a05SBarry Smith else PetscCheck(z >= s || (p <= 1 && bz != DM_BOUNDARY_PERIODIC),PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local z-width of domain z %" PetscInt_FMT " is smaller than stencil width s %" PetscInt_FMT,z,s); 35447c6ae99SBarry Smith zs = 0; 3558865f1eaSKarl Rupp for (i=0; i<(rank/(m*n)); i++) zs += lz[i]; 35647c6ae99SBarry Smith ye = ys + y; 35747c6ae99SBarry Smith xe = xs + x; 35847c6ae99SBarry Smith ze = zs + z; 35947c6ae99SBarry Smith 36088661749SPeter Brune /* determine ghost region (Xs) and region scattered into (IXs) */ 361d9c9ebe5SPeter Brune if (xs-s > 0) { 362d9c9ebe5SPeter Brune Xs = xs - s; IXs = xs - s; 36388661749SPeter Brune } else { 3648865f1eaSKarl Rupp if (bx) Xs = xs - s; 3658865f1eaSKarl Rupp else Xs = 0; 36688661749SPeter Brune IXs = 0; 36788661749SPeter Brune } 368d9c9ebe5SPeter Brune if (xe+s <= M) { 369d9c9ebe5SPeter Brune Xe = xe + s; IXe = xe + s; 37088661749SPeter Brune } else { 37188661749SPeter Brune if (bx) { 372d9c9ebe5SPeter Brune Xs = xs - s; Xe = xe + s; 3738865f1eaSKarl Rupp } else Xe = M; 37488661749SPeter Brune IXe = M; 37588661749SPeter Brune } 37647c6ae99SBarry Smith 377bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC || bx == DM_BOUNDARY_MIRROR) { 378d9c9ebe5SPeter Brune IXs = xs - s; 379d9c9ebe5SPeter Brune IXe = xe + s; 380d9c9ebe5SPeter Brune Xs = xs - s; 381d9c9ebe5SPeter Brune Xe = xe + s; 38288661749SPeter Brune } 38388661749SPeter Brune 384d9c9ebe5SPeter Brune if (ys-s > 0) { 385d9c9ebe5SPeter Brune Ys = ys - s; IYs = ys - s; 38688661749SPeter Brune } else { 3878865f1eaSKarl Rupp if (by) Ys = ys - s; 3888865f1eaSKarl Rupp else Ys = 0; 38988661749SPeter Brune IYs = 0; 39088661749SPeter Brune } 391d9c9ebe5SPeter Brune if (ye+s <= N) { 392d9c9ebe5SPeter Brune Ye = ye + s; IYe = ye + s; 39388661749SPeter Brune } else { 3948865f1eaSKarl Rupp if (by) Ye = ye + s; 3958865f1eaSKarl Rupp else Ye = N; 39688661749SPeter Brune IYe = N; 39788661749SPeter Brune } 39888661749SPeter Brune 399bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_PERIODIC || by == DM_BOUNDARY_MIRROR) { 400d9c9ebe5SPeter Brune IYs = ys - s; 401d9c9ebe5SPeter Brune IYe = ye + s; 402d9c9ebe5SPeter Brune Ys = ys - s; 403d9c9ebe5SPeter Brune Ye = ye + s; 40488661749SPeter Brune } 40588661749SPeter Brune 406d9c9ebe5SPeter Brune if (zs-s > 0) { 407d9c9ebe5SPeter Brune Zs = zs - s; IZs = zs - s; 40888661749SPeter Brune } else { 4098865f1eaSKarl Rupp if (bz) Zs = zs - s; 4108865f1eaSKarl Rupp else Zs = 0; 41188661749SPeter Brune IZs = 0; 41288661749SPeter Brune } 413d9c9ebe5SPeter Brune if (ze+s <= P) { 414d9c9ebe5SPeter Brune Ze = ze + s; IZe = ze + s; 41588661749SPeter Brune } else { 4168865f1eaSKarl Rupp if (bz) Ze = ze + s; 4178865f1eaSKarl Rupp else Ze = P; 41888661749SPeter Brune IZe = P; 41988661749SPeter Brune } 42088661749SPeter Brune 421bff4a2f0SMatthew G. Knepley if (bz == DM_BOUNDARY_PERIODIC || bz == DM_BOUNDARY_MIRROR) { 422d9c9ebe5SPeter Brune IZs = zs - s; 423d9c9ebe5SPeter Brune IZe = ze + s; 424d9c9ebe5SPeter Brune Zs = zs - s; 425d9c9ebe5SPeter Brune Ze = ze + s; 42688661749SPeter Brune } 42747c6ae99SBarry Smith 42847c6ae99SBarry Smith /* Resize all X parameters to reflect w */ 429d9c9ebe5SPeter Brune s_x = s; 430d9c9ebe5SPeter Brune s_y = s; 431d9c9ebe5SPeter Brune s_z = s; 43247c6ae99SBarry Smith 43347c6ae99SBarry Smith /* determine starting point of each processor */ 43447c6ae99SBarry Smith nn = x*y*z; 4359566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(size+1,&bases,size,&ldims)); 4369566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm)); 43747c6ae99SBarry Smith bases[0] = 0; 4388865f1eaSKarl Rupp for (i=1; i<=size; i++) bases[i] = ldims[i-1]; 4398865f1eaSKarl Rupp for (i=1; i<=size; i++) bases[i] += bases[i-1]; 440ce00eea3SSatish Balay base = bases[rank]*dof; 44147c6ae99SBarry Smith 44247c6ae99SBarry Smith /* allocate the base parallel and sequential vectors */ 443ce00eea3SSatish Balay dd->Nlocal = x*y*z*dof; 4449566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,NULL,&global)); 445ce00eea3SSatish Balay dd->nlocal = (Xe-Xs)*(Ye-Ys)*(Ze-Zs)*dof; 4469566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,NULL,&local)); 44747c6ae99SBarry Smith 4484104a7a0SPatrick Sanan /* generate global to local vector scatter and local to global mapping*/ 44947c6ae99SBarry Smith 450ce00eea3SSatish Balay /* global to local must include ghost points within the domain, 451ce00eea3SSatish Balay but not ghost points outside the domain that aren't periodic */ 4529566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((IXe-IXs)*(IYe-IYs)*(IZe-IZs),&idx)); 453d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_BOX) { 454ce00eea3SSatish Balay left = IXs - Xs; right = left + (IXe-IXs); 455ce00eea3SSatish Balay bottom = IYs - Ys; top = bottom + (IYe-IYs); 456ce00eea3SSatish Balay down = IZs - Zs; up = down + (IZe-IZs); 457ce00eea3SSatish Balay count = 0; 458ce00eea3SSatish Balay for (i=down; i<up; i++) { 459ce00eea3SSatish Balay for (j=bottom; j<top; j++) { 460ce00eea3SSatish Balay for (k=left; k<right; k++) { 461ce00eea3SSatish Balay idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 462ce00eea3SSatish Balay } 463ce00eea3SSatish Balay } 464ce00eea3SSatish Balay } 4659566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to)); 46647c6ae99SBarry Smith } else { 46747c6ae99SBarry Smith /* This is way ugly! We need to list the funny cross type region */ 468ce00eea3SSatish Balay left = xs - Xs; right = left + x; 46947c6ae99SBarry Smith bottom = ys - Ys; top = bottom + y; 47047c6ae99SBarry Smith down = zs - Zs; up = down + z; 47147c6ae99SBarry Smith count = 0; 472a5b23f4aSJose E. Roman /* the bottom chunk */ 473ce00eea3SSatish Balay for (i=(IZs-Zs); i<down; i++) { 47447c6ae99SBarry Smith for (j=bottom; j<top; j++) { 475ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 47647c6ae99SBarry Smith } 47747c6ae99SBarry Smith } 47847c6ae99SBarry Smith /* the middle piece */ 47947c6ae99SBarry Smith for (i=down; i<up; i++) { 48047c6ae99SBarry Smith /* front */ 481ce00eea3SSatish Balay for (j=(IYs-Ys); j<bottom; j++) { 482ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 48347c6ae99SBarry Smith } 48447c6ae99SBarry Smith /* middle */ 48547c6ae99SBarry Smith for (j=bottom; j<top; j++) { 486ce00eea3SSatish Balay for (k=IXs-Xs; k<IXe-Xs; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 48747c6ae99SBarry Smith } 48847c6ae99SBarry Smith /* back */ 489ce00eea3SSatish Balay for (j=top; j<top+IYe-ye; j++) { 490ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 49147c6ae99SBarry Smith } 49247c6ae99SBarry Smith } 49347c6ae99SBarry Smith /* the top piece */ 494ce00eea3SSatish Balay for (i=up; i<up+IZe-ze; i++) { 49547c6ae99SBarry Smith for (j=bottom; j<top; j++) { 496ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 49747c6ae99SBarry Smith } 49847c6ae99SBarry Smith } 4999566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to)); 50047c6ae99SBarry Smith } 50147c6ae99SBarry Smith 50247c6ae99SBarry Smith /* determine who lies on each side of use stored in n24 n25 n26 50347c6ae99SBarry Smith n21 n22 n23 50447c6ae99SBarry Smith n18 n19 n20 50547c6ae99SBarry Smith 50647c6ae99SBarry Smith n15 n16 n17 50747c6ae99SBarry Smith n12 n14 50847c6ae99SBarry Smith n9 n10 n11 50947c6ae99SBarry Smith 51047c6ae99SBarry Smith n6 n7 n8 51147c6ae99SBarry Smith n3 n4 n5 51247c6ae99SBarry Smith n0 n1 n2 51347c6ae99SBarry Smith */ 51447c6ae99SBarry Smith 51547c6ae99SBarry Smith /* Solve for X,Y, and Z Periodic Case First, Then Modify Solution */ 51647c6ae99SBarry Smith /* Assume Nodes are Internal to the Cube */ 51747c6ae99SBarry Smith n0 = rank - m*n - m - 1; 51847c6ae99SBarry Smith n1 = rank - m*n - m; 51947c6ae99SBarry Smith n2 = rank - m*n - m + 1; 52047c6ae99SBarry Smith n3 = rank - m*n -1; 52147c6ae99SBarry Smith n4 = rank - m*n; 52247c6ae99SBarry Smith n5 = rank - m*n + 1; 52347c6ae99SBarry Smith n6 = rank - m*n + m - 1; 52447c6ae99SBarry Smith n7 = rank - m*n + m; 52547c6ae99SBarry Smith n8 = rank - m*n + m + 1; 52647c6ae99SBarry Smith 52747c6ae99SBarry Smith n9 = rank - m - 1; 52847c6ae99SBarry Smith n10 = rank - m; 52947c6ae99SBarry Smith n11 = rank - m + 1; 53047c6ae99SBarry Smith n12 = rank - 1; 53147c6ae99SBarry Smith n14 = rank + 1; 53247c6ae99SBarry Smith n15 = rank + m - 1; 53347c6ae99SBarry Smith n16 = rank + m; 53447c6ae99SBarry Smith n17 = rank + m + 1; 53547c6ae99SBarry Smith 53647c6ae99SBarry Smith n18 = rank + m*n - m - 1; 53747c6ae99SBarry Smith n19 = rank + m*n - m; 53847c6ae99SBarry Smith n20 = rank + m*n - m + 1; 53947c6ae99SBarry Smith n21 = rank + m*n - 1; 54047c6ae99SBarry Smith n22 = rank + m*n; 54147c6ae99SBarry Smith n23 = rank + m*n + 1; 54247c6ae99SBarry Smith n24 = rank + m*n + m - 1; 54347c6ae99SBarry Smith n25 = rank + m*n + m; 54447c6ae99SBarry Smith n26 = rank + m*n + m + 1; 54547c6ae99SBarry Smith 54647c6ae99SBarry Smith /* Assume Pieces are on Faces of Cube */ 54747c6ae99SBarry Smith 54847c6ae99SBarry Smith if (xs == 0) { /* First assume not corner or edge */ 54947c6ae99SBarry Smith n0 = rank -1 - (m*n); 55047c6ae99SBarry Smith n3 = rank + m -1 - (m*n); 55147c6ae99SBarry Smith n6 = rank + 2*m -1 - (m*n); 55247c6ae99SBarry Smith n9 = rank -1; 55347c6ae99SBarry Smith n12 = rank + m -1; 55447c6ae99SBarry Smith n15 = rank + 2*m -1; 55547c6ae99SBarry Smith n18 = rank -1 + (m*n); 55647c6ae99SBarry Smith n21 = rank + m -1 + (m*n); 55747c6ae99SBarry Smith n24 = rank + 2*m -1 + (m*n); 55847c6ae99SBarry Smith } 55947c6ae99SBarry Smith 560ce00eea3SSatish Balay if (xe == M) { /* First assume not corner or edge */ 56147c6ae99SBarry Smith n2 = rank -2*m +1 - (m*n); 56247c6ae99SBarry Smith n5 = rank - m +1 - (m*n); 56347c6ae99SBarry Smith n8 = rank +1 - (m*n); 56447c6ae99SBarry Smith n11 = rank -2*m +1; 56547c6ae99SBarry Smith n14 = rank - m +1; 56647c6ae99SBarry Smith n17 = rank +1; 56747c6ae99SBarry Smith n20 = rank -2*m +1 + (m*n); 56847c6ae99SBarry Smith n23 = rank - m +1 + (m*n); 56947c6ae99SBarry Smith n26 = rank +1 + (m*n); 57047c6ae99SBarry Smith } 57147c6ae99SBarry Smith 57247c6ae99SBarry Smith if (ys==0) { /* First assume not corner or edge */ 57347c6ae99SBarry Smith n0 = rank + m * (n-1) -1 - (m*n); 57447c6ae99SBarry Smith n1 = rank + m * (n-1) - (m*n); 57547c6ae99SBarry Smith n2 = rank + m * (n-1) +1 - (m*n); 57647c6ae99SBarry Smith n9 = rank + m * (n-1) -1; 57747c6ae99SBarry Smith n10 = rank + m * (n-1); 57847c6ae99SBarry Smith n11 = rank + m * (n-1) +1; 57947c6ae99SBarry Smith n18 = rank + m * (n-1) -1 + (m*n); 58047c6ae99SBarry Smith n19 = rank + m * (n-1) + (m*n); 58147c6ae99SBarry Smith n20 = rank + m * (n-1) +1 + (m*n); 58247c6ae99SBarry Smith } 58347c6ae99SBarry Smith 58447c6ae99SBarry Smith if (ye == N) { /* First assume not corner or edge */ 58547c6ae99SBarry Smith n6 = rank - m * (n-1) -1 - (m*n); 58647c6ae99SBarry Smith n7 = rank - m * (n-1) - (m*n); 58747c6ae99SBarry Smith n8 = rank - m * (n-1) +1 - (m*n); 58847c6ae99SBarry Smith n15 = rank - m * (n-1) -1; 58947c6ae99SBarry Smith n16 = rank - m * (n-1); 59047c6ae99SBarry Smith n17 = rank - m * (n-1) +1; 59147c6ae99SBarry Smith n24 = rank - m * (n-1) -1 + (m*n); 59247c6ae99SBarry Smith n25 = rank - m * (n-1) + (m*n); 59347c6ae99SBarry Smith n26 = rank - m * (n-1) +1 + (m*n); 59447c6ae99SBarry Smith } 59547c6ae99SBarry Smith 59647c6ae99SBarry Smith if (zs == 0) { /* First assume not corner or edge */ 59747c6ae99SBarry Smith n0 = size - (m*n) + rank - m - 1; 59847c6ae99SBarry Smith n1 = size - (m*n) + rank - m; 59947c6ae99SBarry Smith n2 = size - (m*n) + rank - m + 1; 60047c6ae99SBarry Smith n3 = size - (m*n) + rank - 1; 60147c6ae99SBarry Smith n4 = size - (m*n) + rank; 60247c6ae99SBarry Smith n5 = size - (m*n) + rank + 1; 60347c6ae99SBarry Smith n6 = size - (m*n) + rank + m - 1; 60447c6ae99SBarry Smith n7 = size - (m*n) + rank + m; 60547c6ae99SBarry Smith n8 = size - (m*n) + rank + m + 1; 60647c6ae99SBarry Smith } 60747c6ae99SBarry Smith 60847c6ae99SBarry Smith if (ze == P) { /* First assume not corner or edge */ 60947c6ae99SBarry Smith n18 = (m*n) - (size-rank) - m - 1; 61047c6ae99SBarry Smith n19 = (m*n) - (size-rank) - m; 61147c6ae99SBarry Smith n20 = (m*n) - (size-rank) - m + 1; 61247c6ae99SBarry Smith n21 = (m*n) - (size-rank) - 1; 61347c6ae99SBarry Smith n22 = (m*n) - (size-rank); 61447c6ae99SBarry Smith n23 = (m*n) - (size-rank) + 1; 61547c6ae99SBarry Smith n24 = (m*n) - (size-rank) + m - 1; 61647c6ae99SBarry Smith n25 = (m*n) - (size-rank) + m; 61747c6ae99SBarry Smith n26 = (m*n) - (size-rank) + m + 1; 61847c6ae99SBarry Smith } 61947c6ae99SBarry Smith 62047c6ae99SBarry Smith if ((xs==0) && (zs==0)) { /* Assume an edge, not corner */ 62147c6ae99SBarry Smith n0 = size - m*n + rank + m-1 - m; 62247c6ae99SBarry Smith n3 = size - m*n + rank + m-1; 62347c6ae99SBarry Smith n6 = size - m*n + rank + m-1 + m; 62447c6ae99SBarry Smith } 62547c6ae99SBarry Smith 62647c6ae99SBarry Smith if ((xs==0) && (ze==P)) { /* Assume an edge, not corner */ 62747c6ae99SBarry Smith n18 = m*n - (size - rank) + m-1 - m; 62847c6ae99SBarry Smith n21 = m*n - (size - rank) + m-1; 62947c6ae99SBarry Smith n24 = m*n - (size - rank) + m-1 + m; 63047c6ae99SBarry Smith } 63147c6ae99SBarry Smith 63247c6ae99SBarry Smith if ((xs==0) && (ys==0)) { /* Assume an edge, not corner */ 63347c6ae99SBarry Smith n0 = rank + m*n -1 - m*n; 63447c6ae99SBarry Smith n9 = rank + m*n -1; 63547c6ae99SBarry Smith n18 = rank + m*n -1 + m*n; 63647c6ae99SBarry Smith } 63747c6ae99SBarry Smith 63847c6ae99SBarry Smith if ((xs==0) && (ye==N)) { /* Assume an edge, not corner */ 63947c6ae99SBarry Smith n6 = rank - m*(n-1) + m-1 - m*n; 64047c6ae99SBarry Smith n15 = rank - m*(n-1) + m-1; 64147c6ae99SBarry Smith n24 = rank - m*(n-1) + m-1 + m*n; 64247c6ae99SBarry Smith } 64347c6ae99SBarry Smith 644ce00eea3SSatish Balay if ((xe==M) && (zs==0)) { /* Assume an edge, not corner */ 64547c6ae99SBarry Smith n2 = size - (m*n-rank) - (m-1) - m; 64647c6ae99SBarry Smith n5 = size - (m*n-rank) - (m-1); 64747c6ae99SBarry Smith n8 = size - (m*n-rank) - (m-1) + m; 64847c6ae99SBarry Smith } 64947c6ae99SBarry Smith 650ce00eea3SSatish Balay if ((xe==M) && (ze==P)) { /* Assume an edge, not corner */ 65147c6ae99SBarry Smith n20 = m*n - (size - rank) - (m-1) - m; 65247c6ae99SBarry Smith n23 = m*n - (size - rank) - (m-1); 65347c6ae99SBarry Smith n26 = m*n - (size - rank) - (m-1) + m; 65447c6ae99SBarry Smith } 65547c6ae99SBarry Smith 656ce00eea3SSatish Balay if ((xe==M) && (ys==0)) { /* Assume an edge, not corner */ 65747c6ae99SBarry Smith n2 = rank + m*(n-1) - (m-1) - m*n; 65847c6ae99SBarry Smith n11 = rank + m*(n-1) - (m-1); 65947c6ae99SBarry Smith n20 = rank + m*(n-1) - (m-1) + m*n; 66047c6ae99SBarry Smith } 66147c6ae99SBarry Smith 662ce00eea3SSatish Balay if ((xe==M) && (ye==N)) { /* Assume an edge, not corner */ 66347c6ae99SBarry Smith n8 = rank - m*n +1 - m*n; 66447c6ae99SBarry Smith n17 = rank - m*n +1; 66547c6ae99SBarry Smith n26 = rank - m*n +1 + m*n; 66647c6ae99SBarry Smith } 66747c6ae99SBarry Smith 66847c6ae99SBarry Smith if ((ys==0) && (zs==0)) { /* Assume an edge, not corner */ 66947c6ae99SBarry Smith n0 = size - m + rank -1; 67047c6ae99SBarry Smith n1 = size - m + rank; 67147c6ae99SBarry Smith n2 = size - m + rank +1; 67247c6ae99SBarry Smith } 67347c6ae99SBarry Smith 67447c6ae99SBarry Smith if ((ys==0) && (ze==P)) { /* Assume an edge, not corner */ 67547c6ae99SBarry Smith n18 = m*n - (size - rank) + m*(n-1) -1; 67647c6ae99SBarry Smith n19 = m*n - (size - rank) + m*(n-1); 67747c6ae99SBarry Smith n20 = m*n - (size - rank) + m*(n-1) +1; 67847c6ae99SBarry Smith } 67947c6ae99SBarry Smith 68047c6ae99SBarry Smith if ((ye==N) && (zs==0)) { /* Assume an edge, not corner */ 68147c6ae99SBarry Smith n6 = size - (m*n-rank) - m * (n-1) -1; 68247c6ae99SBarry Smith n7 = size - (m*n-rank) - m * (n-1); 68347c6ae99SBarry Smith n8 = size - (m*n-rank) - m * (n-1) +1; 68447c6ae99SBarry Smith } 68547c6ae99SBarry Smith 68647c6ae99SBarry Smith if ((ye==N) && (ze==P)) { /* Assume an edge, not corner */ 68747c6ae99SBarry Smith n24 = rank - (size-m) -1; 68847c6ae99SBarry Smith n25 = rank - (size-m); 68947c6ae99SBarry Smith n26 = rank - (size-m) +1; 69047c6ae99SBarry Smith } 69147c6ae99SBarry Smith 69247c6ae99SBarry Smith /* Check for Corners */ 6938865f1eaSKarl Rupp if ((xs==0) && (ys==0) && (zs==0)) n0 = size -1; 6948865f1eaSKarl Rupp if ((xs==0) && (ys==0) && (ze==P)) n18 = m*n-1; 6958865f1eaSKarl Rupp if ((xs==0) && (ye==N) && (zs==0)) n6 = (size-1)-m*(n-1); 6968865f1eaSKarl Rupp if ((xs==0) && (ye==N) && (ze==P)) n24 = m-1; 6978865f1eaSKarl Rupp if ((xe==M) && (ys==0) && (zs==0)) n2 = size-m; 6988865f1eaSKarl Rupp if ((xe==M) && (ys==0) && (ze==P)) n20 = m*n-m; 6998865f1eaSKarl Rupp if ((xe==M) && (ye==N) && (zs==0)) n8 = size-m*n; 7008865f1eaSKarl Rupp if ((xe==M) && (ye==N) && (ze==P)) n26 = 0; 70147c6ae99SBarry Smith 70247c6ae99SBarry Smith /* Check for when not X,Y, and Z Periodic */ 70347c6ae99SBarry Smith 70447c6ae99SBarry Smith /* If not X periodic */ 705bff4a2f0SMatthew G. Knepley if (bx != DM_BOUNDARY_PERIODIC) { 7068865f1eaSKarl Rupp if (xs==0) n0 = n3 = n6 = n9 = n12 = n15 = n18 = n21 = n24 = -2; 7078865f1eaSKarl Rupp if (xe==M) n2 = n5 = n8 = n11 = n14 = n17 = n20 = n23 = n26 = -2; 70847c6ae99SBarry Smith } 70947c6ae99SBarry Smith 71047c6ae99SBarry Smith /* If not Y periodic */ 711bff4a2f0SMatthew G. Knepley if (by != DM_BOUNDARY_PERIODIC) { 7128865f1eaSKarl Rupp if (ys==0) n0 = n1 = n2 = n9 = n10 = n11 = n18 = n19 = n20 = -2; 7138865f1eaSKarl Rupp if (ye==N) n6 = n7 = n8 = n15 = n16 = n17 = n24 = n25 = n26 = -2; 71447c6ae99SBarry Smith } 71547c6ae99SBarry Smith 71647c6ae99SBarry Smith /* If not Z periodic */ 717bff4a2f0SMatthew G. Knepley if (bz != DM_BOUNDARY_PERIODIC) { 7188865f1eaSKarl Rupp if (zs==0) n0 = n1 = n2 = n3 = n4 = n5 = n6 = n7 = n8 = -2; 7198865f1eaSKarl Rupp if (ze==P) n18 = n19 = n20 = n21 = n22 = n23 = n24 = n25 = n26 = -2; 72047c6ae99SBarry Smith } 72147c6ae99SBarry Smith 7229566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(27,&dd->neighbors)); 7238865f1eaSKarl Rupp 72447c6ae99SBarry Smith dd->neighbors[0] = n0; 72547c6ae99SBarry Smith dd->neighbors[1] = n1; 72647c6ae99SBarry Smith dd->neighbors[2] = n2; 72747c6ae99SBarry Smith dd->neighbors[3] = n3; 72847c6ae99SBarry Smith dd->neighbors[4] = n4; 72947c6ae99SBarry Smith dd->neighbors[5] = n5; 73047c6ae99SBarry Smith dd->neighbors[6] = n6; 73147c6ae99SBarry Smith dd->neighbors[7] = n7; 73247c6ae99SBarry Smith dd->neighbors[8] = n8; 73347c6ae99SBarry Smith dd->neighbors[9] = n9; 73447c6ae99SBarry Smith dd->neighbors[10] = n10; 73547c6ae99SBarry Smith dd->neighbors[11] = n11; 73647c6ae99SBarry Smith dd->neighbors[12] = n12; 73747c6ae99SBarry Smith dd->neighbors[13] = rank; 73847c6ae99SBarry Smith dd->neighbors[14] = n14; 73947c6ae99SBarry Smith dd->neighbors[15] = n15; 74047c6ae99SBarry Smith dd->neighbors[16] = n16; 74147c6ae99SBarry Smith dd->neighbors[17] = n17; 74247c6ae99SBarry Smith dd->neighbors[18] = n18; 74347c6ae99SBarry Smith dd->neighbors[19] = n19; 74447c6ae99SBarry Smith dd->neighbors[20] = n20; 74547c6ae99SBarry Smith dd->neighbors[21] = n21; 74647c6ae99SBarry Smith dd->neighbors[22] = n22; 74747c6ae99SBarry Smith dd->neighbors[23] = n23; 74847c6ae99SBarry Smith dd->neighbors[24] = n24; 74947c6ae99SBarry Smith dd->neighbors[25] = n25; 75047c6ae99SBarry Smith dd->neighbors[26] = n26; 75147c6ae99SBarry Smith 75247c6ae99SBarry Smith /* If star stencil then delete the corner neighbors */ 753d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_STAR) { 75447c6ae99SBarry Smith /* save information about corner neighbors */ 75547c6ae99SBarry Smith sn0 = n0; sn1 = n1; sn2 = n2; sn3 = n3; sn5 = n5; sn6 = n6; sn7 = n7; 75647c6ae99SBarry Smith sn8 = n8; sn9 = n9; sn11 = n11; sn15 = n15; sn17 = n17; sn18 = n18; 75747c6ae99SBarry Smith sn19 = n19; sn20 = n20; sn21 = n21; sn23 = n23; sn24 = n24; sn25 = n25; 75847c6ae99SBarry Smith sn26 = n26; 7598865f1eaSKarl Rupp n0 = n1 = n2 = n3 = n5 = n6 = n7 = n8 = n9 = n11 = n15 = n17 = n18 = n19 = n20 = n21 = n23 = n24 = n25 = n26 = -1; 76047c6ae99SBarry Smith } 76147c6ae99SBarry Smith 7629566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((Xe-Xs)*(Ye-Ys)*(Ze-Zs),&idx)); 76347c6ae99SBarry Smith 76447c6ae99SBarry Smith nn = 0; 76547c6ae99SBarry Smith /* Bottom Level */ 76647c6ae99SBarry Smith for (k=0; k<s_z; k++) { 76747c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 76847c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 769ce00eea3SSatish Balay x_t = lx[n0 % m]; 77047c6ae99SBarry Smith y_t = ly[(n0 % (m*n))/m]; 77147c6ae99SBarry Smith z_t = lz[n0 / (m*n)]; 77247c6ae99SBarry 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; 7738865f1eaSKarl 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 */ 7748865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 77547c6ae99SBarry Smith } 77647c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 77747c6ae99SBarry Smith x_t = x; 77847c6ae99SBarry Smith y_t = ly[(n1 % (m*n))/m]; 77947c6ae99SBarry Smith z_t = lz[n1 / (m*n)]; 78047c6ae99SBarry 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; 7818865f1eaSKarl Rupp if (twod && (s_t < 0)) s_t = bases[n1] + x_t*y_t*z_t - (s_y+1-i)*x_t; /* 2D case */ 7828865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 78347c6ae99SBarry Smith } 78447c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 785ce00eea3SSatish Balay x_t = lx[n2 % m]; 78647c6ae99SBarry Smith y_t = ly[(n2 % (m*n))/m]; 78747c6ae99SBarry Smith z_t = lz[n2 / (m*n)]; 78847c6ae99SBarry 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; 7898865f1eaSKarl Rupp if (twod && (s_t < 0)) s_t = bases[n2] + x_t*y_t*z_t - (s_y+1-i)*x_t; /* 2D case */ 7908865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 79147c6ae99SBarry Smith } 79247c6ae99SBarry Smith } 79347c6ae99SBarry Smith 79447c6ae99SBarry Smith for (i=0; i<y; i++) { 79547c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 796ce00eea3SSatish Balay x_t = lx[n3 % m]; 79747c6ae99SBarry Smith y_t = y; 79847c6ae99SBarry Smith z_t = lz[n3 / (m*n)]; 79947c6ae99SBarry Smith s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8008865f1eaSKarl 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 */ 8018865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 80247c6ae99SBarry Smith } 80347c6ae99SBarry Smith 80447c6ae99SBarry Smith if (n4 >= 0) { /* middle */ 80547c6ae99SBarry Smith x_t = x; 80647c6ae99SBarry Smith y_t = y; 80747c6ae99SBarry Smith z_t = lz[n4 / (m*n)]; 80847c6ae99SBarry Smith s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8098865f1eaSKarl 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 */ 8108865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 811bff4a2f0SMatthew G. Knepley } else if (bz == DM_BOUNDARY_MIRROR) { 812a6fd6b0aSBarry Smith for (j=0; j<x; j++) idx[nn++] = bases[rank] + j + x*i + (s_z - k - 1)*x*y; 81347c6ae99SBarry Smith } 81447c6ae99SBarry Smith 81547c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 816ce00eea3SSatish Balay x_t = lx[n5 % m]; 81747c6ae99SBarry Smith y_t = y; 81847c6ae99SBarry Smith z_t = lz[n5 / (m*n)]; 81947c6ae99SBarry Smith s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8208865f1eaSKarl 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 */ 8218865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 82247c6ae99SBarry Smith } 82347c6ae99SBarry Smith } 82447c6ae99SBarry Smith 82547c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 82647c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 827ce00eea3SSatish Balay x_t = lx[n6 % m]; 82847c6ae99SBarry Smith y_t = ly[(n6 % (m*n))/m]; 82947c6ae99SBarry Smith z_t = lz[n6 / (m*n)]; 83047c6ae99SBarry Smith s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8318865f1eaSKarl 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 */ 8328865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 83347c6ae99SBarry Smith } 83447c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 83547c6ae99SBarry Smith x_t = x; 83647c6ae99SBarry Smith y_t = ly[(n7 % (m*n))/m]; 83747c6ae99SBarry Smith z_t = lz[n7 / (m*n)]; 83847c6ae99SBarry Smith s_t = bases[n7] + (i-1)*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[n7] + (i-1)*x_t + x_t*y_t*z_t - x_t*y_t; /* 2D case */ 8408865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 84147c6ae99SBarry Smith } 84247c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 843ce00eea3SSatish Balay x_t = lx[n8 % m]; 84447c6ae99SBarry Smith y_t = ly[(n8 % (m*n))/m]; 84547c6ae99SBarry Smith z_t = lz[n8 / (m*n)]; 84647c6ae99SBarry Smith s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8478865f1eaSKarl 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 */ 8488865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 84947c6ae99SBarry Smith } 85047c6ae99SBarry Smith } 85147c6ae99SBarry Smith } 85247c6ae99SBarry Smith 85347c6ae99SBarry Smith /* Middle Level */ 85447c6ae99SBarry Smith for (k=0; k<z; k++) { 85547c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 85647c6ae99SBarry Smith if (n9 >= 0) { /* left below */ 857ce00eea3SSatish Balay x_t = lx[n9 % m]; 85847c6ae99SBarry Smith y_t = ly[(n9 % (m*n))/m]; 85947c6ae99SBarry Smith /* z_t = z; */ 86047c6ae99SBarry Smith s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 8618865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 86247c6ae99SBarry Smith } 86347c6ae99SBarry Smith if (n10 >= 0) { /* directly below */ 86447c6ae99SBarry Smith x_t = x; 86547c6ae99SBarry Smith y_t = ly[(n10 % (m*n))/m]; 86647c6ae99SBarry Smith /* z_t = z; */ 86747c6ae99SBarry Smith s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 8688865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 869bff4a2f0SMatthew G. Knepley } else if (by == DM_BOUNDARY_MIRROR) { 870a6fd6b0aSBarry Smith for (j=0; j<x; j++) idx[nn++] = bases[rank] + j + k*x*y + (s_y - i)*x; 87147c6ae99SBarry Smith } 87247c6ae99SBarry Smith if (n11 >= 0) { /* right below */ 873ce00eea3SSatish Balay x_t = lx[n11 % m]; 87447c6ae99SBarry Smith y_t = ly[(n11 % (m*n))/m]; 87547c6ae99SBarry Smith /* z_t = z; */ 87647c6ae99SBarry Smith s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 8778865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 87847c6ae99SBarry Smith } 87947c6ae99SBarry Smith } 88047c6ae99SBarry Smith 88147c6ae99SBarry Smith for (i=0; i<y; i++) { 88247c6ae99SBarry Smith if (n12 >= 0) { /* directly left */ 883ce00eea3SSatish Balay x_t = lx[n12 % m]; 88447c6ae99SBarry Smith y_t = y; 88547c6ae99SBarry Smith /* z_t = z; */ 88647c6ae99SBarry Smith s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t; 8878865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 888bff4a2f0SMatthew G. Knepley } else if (bx == DM_BOUNDARY_MIRROR) { 889a6fd6b0aSBarry Smith for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + s_x - j - 1 + k*x*y + i*x; 89047c6ae99SBarry Smith } 89147c6ae99SBarry Smith 89247c6ae99SBarry Smith /* Interior */ 89347c6ae99SBarry Smith s_t = bases[rank] + i*x + k*x*y; 8948865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = s_t++; 89547c6ae99SBarry Smith 89647c6ae99SBarry Smith if (n14 >= 0) { /* directly right */ 897ce00eea3SSatish Balay x_t = lx[n14 % m]; 89847c6ae99SBarry Smith y_t = y; 89947c6ae99SBarry Smith /* z_t = z; */ 90047c6ae99SBarry Smith s_t = bases[n14] + i*x_t + k*x_t*y_t; 9018865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 902bff4a2f0SMatthew G. Knepley } else if (bx == DM_BOUNDARY_MIRROR) { 903a6fd6b0aSBarry Smith for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x - j - 1 + k*x*y + i*x; 90447c6ae99SBarry Smith } 90547c6ae99SBarry Smith } 90647c6ae99SBarry Smith 90747c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 90847c6ae99SBarry Smith if (n15 >= 0) { /* left above */ 909ce00eea3SSatish Balay x_t = lx[n15 % m]; 91047c6ae99SBarry Smith y_t = ly[(n15 % (m*n))/m]; 91147c6ae99SBarry Smith /* z_t = z; */ 91247c6ae99SBarry Smith s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t; 9138865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 91447c6ae99SBarry Smith } 91547c6ae99SBarry Smith if (n16 >= 0) { /* directly above */ 91647c6ae99SBarry Smith x_t = x; 91747c6ae99SBarry Smith y_t = ly[(n16 % (m*n))/m]; 91847c6ae99SBarry Smith /* z_t = z; */ 91947c6ae99SBarry Smith s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t; 9208865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 921bff4a2f0SMatthew G. Knepley } else if (by == DM_BOUNDARY_MIRROR) { 922a6fd6b0aSBarry Smith for (j=0; j<x; j++) idx[nn++] = bases[rank] + j + k*x*y + (y-i)*x; 92347c6ae99SBarry Smith } 92447c6ae99SBarry Smith if (n17 >= 0) { /* right above */ 925ce00eea3SSatish Balay x_t = lx[n17 % m]; 92647c6ae99SBarry Smith y_t = ly[(n17 % (m*n))/m]; 92747c6ae99SBarry Smith /* z_t = z; */ 92847c6ae99SBarry Smith s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t; 9298865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 93047c6ae99SBarry Smith } 93147c6ae99SBarry Smith } 93247c6ae99SBarry Smith } 93347c6ae99SBarry Smith 93447c6ae99SBarry Smith /* Upper Level */ 93547c6ae99SBarry Smith for (k=0; k<s_z; k++) { 93647c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 93747c6ae99SBarry Smith if (n18 >= 0) { /* left below */ 938ce00eea3SSatish Balay x_t = lx[n18 % m]; 93947c6ae99SBarry Smith y_t = ly[(n18 % (m*n))/m]; 94047c6ae99SBarry Smith /* z_t = lz[n18 / (m*n)]; */ 94147c6ae99SBarry Smith s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 9428865f1eaSKarl 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 */ 9438865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 94447c6ae99SBarry Smith } 94547c6ae99SBarry Smith if (n19 >= 0) { /* directly below */ 94647c6ae99SBarry Smith x_t = x; 94747c6ae99SBarry Smith y_t = ly[(n19 % (m*n))/m]; 94847c6ae99SBarry Smith /* z_t = lz[n19 / (m*n)]; */ 94947c6ae99SBarry Smith s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 9508865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n19] - (s_y+1-i)*x_t + x_t*y_t; /* 2d case */ 9518865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 95247c6ae99SBarry Smith } 95347c6ae99SBarry Smith if (n20 >= 0) { /* right below */ 954ce00eea3SSatish Balay x_t = lx[n20 % m]; 95547c6ae99SBarry Smith y_t = ly[(n20 % (m*n))/m]; 95647c6ae99SBarry Smith /* z_t = lz[n20 / (m*n)]; */ 95747c6ae99SBarry Smith s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 9588865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n20] - (s_y+1-i)*x_t + x_t*y_t; /* 2d case */ 9598865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 96047c6ae99SBarry Smith } 96147c6ae99SBarry Smith } 96247c6ae99SBarry Smith 96347c6ae99SBarry Smith for (i=0; i<y; i++) { 96447c6ae99SBarry Smith if (n21 >= 0) { /* directly left */ 965ce00eea3SSatish Balay x_t = lx[n21 % m]; 96647c6ae99SBarry Smith y_t = y; 96747c6ae99SBarry Smith /* z_t = lz[n21 / (m*n)]; */ 96847c6ae99SBarry Smith s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t; 9698865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n21] + (i+1)*x_t - s_x; /* 2d case */ 9708865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 97147c6ae99SBarry Smith } 97247c6ae99SBarry Smith 97347c6ae99SBarry Smith if (n22 >= 0) { /* middle */ 97447c6ae99SBarry Smith x_t = x; 97547c6ae99SBarry Smith y_t = y; 97647c6ae99SBarry Smith /* z_t = lz[n22 / (m*n)]; */ 97747c6ae99SBarry Smith s_t = bases[n22] + i*x_t + k*x_t*y_t; 9788865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n22] + i*x_t; /* 2d case */ 9798865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 980bff4a2f0SMatthew G. Knepley } else if (bz == DM_BOUNDARY_MIRROR) { 981a6fd6b0aSBarry Smith for (j=0; j<x; j++) idx[nn++] = bases[rank] + j + (z-k-1)*x*y + i*x; 98247c6ae99SBarry Smith } 98347c6ae99SBarry Smith 98447c6ae99SBarry Smith if (n23 >= 0) { /* directly right */ 985ce00eea3SSatish Balay x_t = lx[n23 % m]; 98647c6ae99SBarry Smith y_t = y; 98747c6ae99SBarry Smith /* z_t = lz[n23 / (m*n)]; */ 98847c6ae99SBarry Smith s_t = bases[n23] + i*x_t + k*x_t*y_t; 9898865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n23] + i*x_t; /* 2d case */ 9908865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 99147c6ae99SBarry Smith } 99247c6ae99SBarry Smith } 99347c6ae99SBarry Smith 99447c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 99547c6ae99SBarry Smith if (n24 >= 0) { /* left above */ 996ce00eea3SSatish Balay x_t = lx[n24 % m]; 99747c6ae99SBarry Smith y_t = ly[(n24 % (m*n))/m]; 99847c6ae99SBarry Smith /* z_t = lz[n24 / (m*n)]; */ 99947c6ae99SBarry Smith s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t; 10008865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n24] + i*x_t - s_x; /* 2d case */ 10018865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 100247c6ae99SBarry Smith } 100347c6ae99SBarry Smith if (n25 >= 0) { /* directly above */ 100447c6ae99SBarry Smith x_t = x; 100547c6ae99SBarry Smith y_t = ly[(n25 % (m*n))/m]; 100647c6ae99SBarry Smith /* z_t = lz[n25 / (m*n)]; */ 100747c6ae99SBarry Smith s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t; 10088865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n25] + (i-1)*x_t; /* 2d case */ 10098865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 101047c6ae99SBarry Smith } 101147c6ae99SBarry Smith if (n26 >= 0) { /* right above */ 1012ce00eea3SSatish Balay x_t = lx[n26 % m]; 101347c6ae99SBarry Smith y_t = ly[(n26 % (m*n))/m]; 101447c6ae99SBarry Smith /* z_t = lz[n26 / (m*n)]; */ 101547c6ae99SBarry Smith s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t; 10168865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n26] + (i-1)*x_t; /* 2d case */ 10178865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 101847c6ae99SBarry Smith } 101947c6ae99SBarry Smith } 102047c6ae99SBarry Smith } 1021ce00eea3SSatish Balay 10229566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(comm,dof,nn,idx,PETSC_USE_POINTER,&from)); 10239566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(global,from,local,to,>ol)); 10249566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)da,(PetscObject)gtol)); 10259566063dSJacob Faibussowitsch PetscCall(ISDestroy(&to)); 10269566063dSJacob Faibussowitsch PetscCall(ISDestroy(&from)); 102747c6ae99SBarry Smith 1028d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_STAR) { 102947c6ae99SBarry Smith n0 = sn0; n1 = sn1; n2 = sn2; n3 = sn3; n5 = sn5; n6 = sn6; n7 = sn7; 103047c6ae99SBarry Smith n8 = sn8; n9 = sn9; n11 = sn11; n15 = sn15; n17 = sn17; n18 = sn18; 103147c6ae99SBarry Smith n19 = sn19; n20 = sn20; n21 = sn21; n23 = sn23; n24 = sn24; n25 = sn25; 103247c6ae99SBarry Smith n26 = sn26; 1033ce00eea3SSatish Balay } 103447c6ae99SBarry Smith 1035288e7d53SBarry Smith if (((stencil_type == DMDA_STENCIL_STAR) || (bx != DM_BOUNDARY_PERIODIC && bx) || (by != DM_BOUNDARY_PERIODIC && by) || (bz != DM_BOUNDARY_PERIODIC && bz))) { 1036ce00eea3SSatish Balay /* 1037ce00eea3SSatish Balay Recompute the local to global mappings, this time keeping the 1038ce00eea3SSatish Balay information about the cross corner processor numbers. 1039ce00eea3SSatish Balay */ 104047c6ae99SBarry Smith nn = 0; 104147c6ae99SBarry Smith /* Bottom Level */ 104247c6ae99SBarry Smith for (k=0; k<s_z; k++) { 104347c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 104447c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 1045ce00eea3SSatish Balay x_t = lx[n0 % m]; 104647c6ae99SBarry Smith y_t = ly[(n0 % (m*n))/m]; 104747c6ae99SBarry Smith z_t = lz[n0 / (m*n)]; 104847c6ae99SBarry 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; 10498865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1050ce00eea3SSatish Balay } else if (Xs-xs < 0 && Ys-ys < 0 && Zs-zs < 0) { 10518865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 105247c6ae99SBarry Smith } 105347c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 105447c6ae99SBarry Smith x_t = x; 105547c6ae99SBarry Smith y_t = ly[(n1 % (m*n))/m]; 105647c6ae99SBarry Smith z_t = lz[n1 / (m*n)]; 105747c6ae99SBarry 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; 10588865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1059ce00eea3SSatish Balay } else if (Ys-ys < 0 && Zs-zs < 0) { 10608865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 106147c6ae99SBarry Smith } 106247c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 1063ce00eea3SSatish Balay x_t = lx[n2 % m]; 106447c6ae99SBarry Smith y_t = ly[(n2 % (m*n))/m]; 106547c6ae99SBarry Smith z_t = lz[n2 / (m*n)]; 106647c6ae99SBarry 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; 10678865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1068ce00eea3SSatish Balay } else if (xe-Xe < 0 && Ys-ys < 0 && Zs-zs < 0) { 10698865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 107047c6ae99SBarry Smith } 107147c6ae99SBarry Smith } 107247c6ae99SBarry Smith 107347c6ae99SBarry Smith for (i=0; i<y; i++) { 107447c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 1075ce00eea3SSatish Balay x_t = lx[n3 % m]; 107647c6ae99SBarry Smith y_t = y; 107747c6ae99SBarry Smith z_t = lz[n3 / (m*n)]; 107847c6ae99SBarry Smith s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 10798865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1080ce00eea3SSatish Balay } else if (Xs-xs < 0 && Zs-zs < 0) { 10818865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 108247c6ae99SBarry Smith } 108347c6ae99SBarry Smith 108447c6ae99SBarry Smith if (n4 >= 0) { /* middle */ 108547c6ae99SBarry Smith x_t = x; 108647c6ae99SBarry Smith y_t = y; 108747c6ae99SBarry Smith z_t = lz[n4 / (m*n)]; 108847c6ae99SBarry Smith s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 10898865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1090ce00eea3SSatish Balay } else if (Zs-zs < 0) { 1091bff4a2f0SMatthew G. Knepley if (bz == DM_BOUNDARY_MIRROR) { 1092a6fd6b0aSBarry Smith for (j=0; j<x; j++) idx[nn++] = bases[rank] + j + x*i + (s_z - k - 1)*x*y; 1093c2859e5eSBarry Smith } else { 10948865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 109547c6ae99SBarry Smith } 1096c2859e5eSBarry Smith } 109747c6ae99SBarry Smith 109847c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 1099ce00eea3SSatish Balay x_t = lx[n5 % m]; 110047c6ae99SBarry Smith y_t = y; 110147c6ae99SBarry Smith z_t = lz[n5 / (m*n)]; 110247c6ae99SBarry Smith s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 11038865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1104ce00eea3SSatish Balay } else if (xe-Xe < 0 && Zs-zs < 0) { 11058865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 110647c6ae99SBarry Smith } 110747c6ae99SBarry Smith } 110847c6ae99SBarry Smith 110947c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 111047c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 1111ce00eea3SSatish Balay x_t = lx[n6 % m]; 111247c6ae99SBarry Smith y_t = ly[(n6 % (m*n))/m]; 111347c6ae99SBarry Smith z_t = lz[n6 / (m*n)]; 111447c6ae99SBarry Smith s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 11158865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1116ce00eea3SSatish Balay } else if (Xs-xs < 0 && ye-Ye < 0 && Zs-zs < 0) { 11178865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 111847c6ae99SBarry Smith } 111947c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 112047c6ae99SBarry Smith x_t = x; 112147c6ae99SBarry Smith y_t = ly[(n7 % (m*n))/m]; 112247c6ae99SBarry Smith z_t = lz[n7 / (m*n)]; 112347c6ae99SBarry Smith s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 11248865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1125ce00eea3SSatish Balay } else if (ye-Ye < 0 && Zs-zs < 0) { 11268865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 112747c6ae99SBarry Smith } 112847c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 1129ce00eea3SSatish Balay x_t = lx[n8 % m]; 113047c6ae99SBarry Smith y_t = ly[(n8 % (m*n))/m]; 113147c6ae99SBarry Smith z_t = lz[n8 / (m*n)]; 113247c6ae99SBarry Smith s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 11338865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1134ce00eea3SSatish Balay } else if (xe-Xe < 0 && ye-Ye < 0 && Zs-zs < 0) { 11358865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 113647c6ae99SBarry Smith } 113747c6ae99SBarry Smith } 113847c6ae99SBarry Smith } 113947c6ae99SBarry Smith 114047c6ae99SBarry Smith /* Middle Level */ 114147c6ae99SBarry Smith for (k=0; k<z; k++) { 114247c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 114347c6ae99SBarry Smith if (n9 >= 0) { /* left below */ 1144ce00eea3SSatish Balay x_t = lx[n9 % m]; 114547c6ae99SBarry Smith y_t = ly[(n9 % (m*n))/m]; 114647c6ae99SBarry Smith /* z_t = z; */ 114747c6ae99SBarry Smith s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 11488865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1149ce00eea3SSatish Balay } else if (Xs-xs < 0 && Ys-ys < 0) { 11508865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 115147c6ae99SBarry Smith } 115247c6ae99SBarry Smith if (n10 >= 0) { /* directly below */ 115347c6ae99SBarry Smith x_t = x; 115447c6ae99SBarry Smith y_t = ly[(n10 % (m*n))/m]; 115547c6ae99SBarry Smith /* z_t = z; */ 115647c6ae99SBarry Smith s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 11578865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1158ce00eea3SSatish Balay } else if (Ys-ys < 0) { 1159bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_MIRROR) { 1160feb237baSPierre Jolivet for (j=0; j<x; j++) idx[nn++] = bases[rank] + j + k*x*y + (s_y - i)*x; 1161c2859e5eSBarry Smith } else { 11628865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 1163c2859e5eSBarry Smith } 116447c6ae99SBarry Smith } 116547c6ae99SBarry Smith if (n11 >= 0) { /* right below */ 1166ce00eea3SSatish Balay x_t = lx[n11 % m]; 116747c6ae99SBarry Smith y_t = ly[(n11 % (m*n))/m]; 116847c6ae99SBarry Smith /* z_t = z; */ 116947c6ae99SBarry Smith s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 11708865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1171ce00eea3SSatish Balay } else if (xe-Xe < 0 && Ys-ys < 0) { 11728865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 117347c6ae99SBarry Smith } 117447c6ae99SBarry Smith } 117547c6ae99SBarry Smith 117647c6ae99SBarry Smith for (i=0; i<y; i++) { 117747c6ae99SBarry Smith if (n12 >= 0) { /* directly left */ 1178ce00eea3SSatish Balay x_t = lx[n12 % m]; 117947c6ae99SBarry Smith y_t = y; 118047c6ae99SBarry Smith /* z_t = z; */ 118147c6ae99SBarry Smith s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t; 11828865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1183ce00eea3SSatish Balay } else if (Xs-xs < 0) { 1184bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_MIRROR) { 1185a6fd6b0aSBarry Smith for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + s_x - j - 1 + k*x*y + i*x; 1186c2859e5eSBarry Smith } else { 11878865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 118847c6ae99SBarry Smith } 1189c2859e5eSBarry Smith } 119047c6ae99SBarry Smith 119147c6ae99SBarry Smith /* Interior */ 119247c6ae99SBarry Smith s_t = bases[rank] + i*x + k*x*y; 11938865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = s_t++; 119447c6ae99SBarry Smith 119547c6ae99SBarry Smith if (n14 >= 0) { /* directly right */ 1196ce00eea3SSatish Balay x_t = lx[n14 % m]; 119747c6ae99SBarry Smith y_t = y; 119847c6ae99SBarry Smith /* z_t = z; */ 119947c6ae99SBarry Smith s_t = bases[n14] + i*x_t + k*x_t*y_t; 12008865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1201ce00eea3SSatish Balay } else if (xe-Xe < 0) { 1202bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_MIRROR) { 1203a6fd6b0aSBarry Smith for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x - j - 1 + k*x*y + i*x; 1204c2859e5eSBarry Smith } else { 12058865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 120647c6ae99SBarry Smith } 120747c6ae99SBarry Smith } 1208c2859e5eSBarry Smith } 120947c6ae99SBarry Smith 121047c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 121147c6ae99SBarry Smith if (n15 >= 0) { /* left above */ 1212ce00eea3SSatish Balay x_t = lx[n15 % m]; 121347c6ae99SBarry Smith y_t = ly[(n15 % (m*n))/m]; 121447c6ae99SBarry Smith /* z_t = z; */ 121547c6ae99SBarry Smith s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t; 12168865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1217ce00eea3SSatish Balay } else if (Xs-xs < 0 && ye-Ye < 0) { 12188865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 121947c6ae99SBarry Smith } 122047c6ae99SBarry Smith if (n16 >= 0) { /* directly above */ 122147c6ae99SBarry Smith x_t = x; 122247c6ae99SBarry Smith y_t = ly[(n16 % (m*n))/m]; 122347c6ae99SBarry Smith /* z_t = z; */ 122447c6ae99SBarry Smith s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t; 12258865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1226ce00eea3SSatish Balay } else if (ye-Ye < 0) { 1227bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_MIRROR) { 1228a6fd6b0aSBarry Smith for (j=0; j<x; j++) idx[nn++] = bases[rank] + j + k*x*y + (y-i)*x; 1229c2859e5eSBarry Smith } else { 12308865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 123147c6ae99SBarry Smith } 1232c2859e5eSBarry Smith } 123347c6ae99SBarry Smith if (n17 >= 0) { /* right above */ 1234ce00eea3SSatish Balay x_t = lx[n17 % m]; 123547c6ae99SBarry Smith y_t = ly[(n17 % (m*n))/m]; 123647c6ae99SBarry Smith /* z_t = z; */ 123747c6ae99SBarry Smith s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t; 12388865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1239ce00eea3SSatish Balay } else if (xe-Xe < 0 && ye-Ye < 0) { 12408865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 124147c6ae99SBarry Smith } 124247c6ae99SBarry Smith } 124347c6ae99SBarry Smith } 124447c6ae99SBarry Smith 124547c6ae99SBarry Smith /* Upper Level */ 124647c6ae99SBarry Smith for (k=0; k<s_z; k++) { 124747c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 124847c6ae99SBarry Smith if (n18 >= 0) { /* left below */ 1249ce00eea3SSatish Balay x_t = lx[n18 % m]; 125047c6ae99SBarry Smith y_t = ly[(n18 % (m*n))/m]; 125147c6ae99SBarry Smith /* z_t = lz[n18 / (m*n)]; */ 125247c6ae99SBarry Smith s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 12538865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1254ce00eea3SSatish Balay } else if (Xs-xs < 0 && Ys-ys < 0 && ze-Ze < 0) { 12558865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 125647c6ae99SBarry Smith } 125747c6ae99SBarry Smith if (n19 >= 0) { /* directly below */ 125847c6ae99SBarry Smith x_t = x; 125947c6ae99SBarry Smith y_t = ly[(n19 % (m*n))/m]; 126047c6ae99SBarry Smith /* z_t = lz[n19 / (m*n)]; */ 126147c6ae99SBarry Smith s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 12628865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1263ce00eea3SSatish Balay } else if (Ys-ys < 0 && ze-Ze < 0) { 12648865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 126547c6ae99SBarry Smith } 126647c6ae99SBarry Smith if (n20 >= 0) { /* right below */ 1267ce00eea3SSatish Balay x_t = lx[n20 % m]; 126847c6ae99SBarry Smith y_t = ly[(n20 % (m*n))/m]; 126947c6ae99SBarry Smith /* z_t = lz[n20 / (m*n)]; */ 127047c6ae99SBarry Smith s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 12718865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1272ce00eea3SSatish Balay } else if (xe-Xe < 0 && Ys-ys < 0 && ze-Ze < 0) { 12738865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 127447c6ae99SBarry Smith } 127547c6ae99SBarry Smith } 127647c6ae99SBarry Smith 127747c6ae99SBarry Smith for (i=0; i<y; i++) { 127847c6ae99SBarry Smith if (n21 >= 0) { /* directly left */ 1279ce00eea3SSatish Balay x_t = lx[n21 % m]; 128047c6ae99SBarry Smith y_t = y; 128147c6ae99SBarry Smith /* z_t = lz[n21 / (m*n)]; */ 128247c6ae99SBarry Smith s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t; 12838865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1284ce00eea3SSatish Balay } else if (Xs-xs < 0 && ze-Ze < 0) { 12858865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 128647c6ae99SBarry Smith } 128747c6ae99SBarry Smith 128847c6ae99SBarry Smith if (n22 >= 0) { /* middle */ 128947c6ae99SBarry Smith x_t = x; 129047c6ae99SBarry Smith y_t = y; 129147c6ae99SBarry Smith /* z_t = lz[n22 / (m*n)]; */ 129247c6ae99SBarry Smith s_t = bases[n22] + i*x_t + k*x_t*y_t; 12938865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1294ce00eea3SSatish Balay } else if (ze-Ze < 0) { 1295bff4a2f0SMatthew G. Knepley if (bz == DM_BOUNDARY_MIRROR) { 1296a6fd6b0aSBarry Smith for (j=0; j<x; j++) idx[nn++] = bases[rank] + j + (z-k-1)*x*y + i*x; 1297c2859e5eSBarry Smith } else { 12988865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 129947c6ae99SBarry Smith } 1300c2859e5eSBarry Smith } 130147c6ae99SBarry Smith 130247c6ae99SBarry Smith if (n23 >= 0) { /* directly right */ 1303ce00eea3SSatish Balay x_t = lx[n23 % m]; 130447c6ae99SBarry Smith y_t = y; 130547c6ae99SBarry Smith /* z_t = lz[n23 / (m*n)]; */ 130647c6ae99SBarry Smith s_t = bases[n23] + i*x_t + k*x_t*y_t; 13078865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1308ce00eea3SSatish Balay } else if (xe-Xe < 0 && ze-Ze < 0) { 13098865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 131047c6ae99SBarry Smith } 131147c6ae99SBarry Smith } 131247c6ae99SBarry Smith 131347c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 131447c6ae99SBarry Smith if (n24 >= 0) { /* left above */ 1315ce00eea3SSatish Balay x_t = lx[n24 % m]; 131647c6ae99SBarry Smith y_t = ly[(n24 % (m*n))/m]; 131747c6ae99SBarry Smith /* z_t = lz[n24 / (m*n)]; */ 131847c6ae99SBarry Smith s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t; 13198865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1320ce00eea3SSatish Balay } else if (Xs-xs < 0 && ye-Ye < 0 && ze-Ze < 0) { 13218865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 132247c6ae99SBarry Smith } 132347c6ae99SBarry Smith if (n25 >= 0) { /* directly above */ 132447c6ae99SBarry Smith x_t = x; 132547c6ae99SBarry Smith y_t = ly[(n25 % (m*n))/m]; 132647c6ae99SBarry Smith /* z_t = lz[n25 / (m*n)]; */ 132747c6ae99SBarry Smith s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t; 13288865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1329ce00eea3SSatish Balay } else if (ye-Ye < 0 && ze-Ze < 0) { 13308865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 133147c6ae99SBarry Smith } 133247c6ae99SBarry Smith if (n26 >= 0) { /* right above */ 1333ce00eea3SSatish Balay x_t = lx[n26 % m]; 133447c6ae99SBarry Smith y_t = ly[(n26 % (m*n))/m]; 133547c6ae99SBarry Smith /* z_t = lz[n26 / (m*n)]; */ 133647c6ae99SBarry Smith s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t; 13378865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1338ce00eea3SSatish Balay } else if (xe-Xe < 0 && ye-Ye < 0 && ze-Ze < 0) { 13398865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 134047c6ae99SBarry Smith } 134147c6ae99SBarry Smith } 134247c6ae99SBarry Smith } 134347c6ae99SBarry Smith } 134447c6ae99SBarry Smith /* 134547c6ae99SBarry Smith Set the local to global ordering in the global vector, this allows use 134647c6ae99SBarry Smith of VecSetValuesLocal(). 134747c6ae99SBarry Smith */ 13489566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(comm,dof,nn,idx,PETSC_OWN_POINTER,&da->ltogmap)); 13499566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)da,(PetscObject)da->ltogmap)); 135047c6ae99SBarry Smith 13519566063dSJacob Faibussowitsch PetscCall(PetscFree2(bases,ldims)); 1352ce00eea3SSatish Balay dd->m = m; dd->n = n; dd->p = p; 1353ce00eea3SSatish Balay /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */ 1354ce00eea3SSatish Balay dd->xs = xs*dof; dd->xe = xe*dof; dd->ys = ys; dd->ye = ye; dd->zs = zs; dd->ze = ze; 1355ce00eea3SSatish Balay dd->Xs = Xs*dof; dd->Xe = Xe*dof; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = Zs; dd->Ze = Ze; 1356ce00eea3SSatish Balay 13579566063dSJacob Faibussowitsch PetscCall(VecDestroy(&local)); 13589566063dSJacob Faibussowitsch PetscCall(VecDestroy(&global)); 1359ce00eea3SSatish Balay 1360ce00eea3SSatish Balay dd->gtol = gtol; 1361ce00eea3SSatish Balay dd->base = base; 1362ce00eea3SSatish Balay da->ops->view = DMView_DA_3d; 13630298fd71SBarry Smith dd->ltol = NULL; 13640298fd71SBarry Smith dd->ao = NULL; 136547c6ae99SBarry Smith PetscFunctionReturn(0); 136647c6ae99SBarry Smith } 1367564755cdSBarry Smith 136847c6ae99SBarry Smith /*@C 1369aa219208SBarry Smith DMDACreate3d - Creates an object that will manage the communication of three-dimensional 137047c6ae99SBarry Smith regular array data that is distributed across some processors. 137147c6ae99SBarry Smith 1372d083f849SBarry Smith Collective 137347c6ae99SBarry Smith 137447c6ae99SBarry Smith Input Parameters: 137547c6ae99SBarry Smith + comm - MPI communicator 13761321219cSEthan Coon . bx,by,bz - type of ghost nodes the array have. 1377bff4a2f0SMatthew G. Knepley Use one of DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_PERIODIC. 1378aa219208SBarry Smith . stencil_type - Type of stencil (DMDA_STENCIL_STAR or DMDA_STENCIL_BOX) 1379897f7067SBarry Smith . M,N,P - global dimension in each direction of the array 138047c6ae99SBarry Smith . m,n,p - corresponding number of processors in each dimension 138147c6ae99SBarry Smith (or PETSC_DECIDE to have calculated) 138247c6ae99SBarry Smith . dof - number of degrees of freedom per node 138310d7c030SMatthew G Knepley . s - stencil width 138410d7c030SMatthew G Knepley - lx, ly, lz - arrays containing the number of nodes in each cell along 13850298fd71SBarry Smith the x, y, and z coordinates, or NULL. If non-null, these 138647c6ae99SBarry Smith must be of length as m,n,p and the corresponding 138747c6ae99SBarry Smith m,n, or p cannot be PETSC_DECIDE. Sum of the lx[] entries must be M, sum of 138847c6ae99SBarry Smith the ly[] must N, sum of the lz[] must be P 138947c6ae99SBarry Smith 139047c6ae99SBarry Smith Output Parameter: 139147c6ae99SBarry Smith . da - the resulting distributed array object 139247c6ae99SBarry Smith 139347c6ae99SBarry Smith Options Database Key: 1394706a11cbSBarry Smith + -dm_view - Calls DMView() at the conclusion of DMDACreate3d() 1395e43dc8daSBarry Smith . -da_grid_x <nx> - number of grid points in x direction 1396e43dc8daSBarry Smith . -da_grid_y <ny> - number of grid points in y direction 1397e43dc8daSBarry Smith . -da_grid_z <nz> - number of grid points in z direction 139847c6ae99SBarry Smith . -da_processors_x <MX> - number of processors in x direction 139947c6ae99SBarry Smith . -da_processors_y <MY> - number of processors in y direction 140047c6ae99SBarry Smith . -da_processors_z <MZ> - number of processors in z direction 1401e0f5d30fSBarry Smith . -da_refine_x <rx> - refinement ratio in x direction 1402e0f5d30fSBarry Smith . -da_refine_y <ry> - refinement ratio in y direction 1403e0f5d30fSBarry Smith . -da_refine_z <rz>- refinement ratio in z directio 1404e43dc8daSBarry Smith - -da_refine <n> - refine the DMDA n times before creating it 140547c6ae99SBarry Smith 140647c6ae99SBarry Smith Level: beginner 140747c6ae99SBarry Smith 140847c6ae99SBarry Smith Notes: 1409aa219208SBarry Smith The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the 1410aa219208SBarry Smith standard 7-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes 141147c6ae99SBarry Smith the standard 27-pt stencil. 141247c6ae99SBarry Smith 1413aa219208SBarry Smith The array data itself is NOT stored in the DMDA, it is stored in Vec objects; 1414564755cdSBarry Smith The appropriate vector objects can be obtained with calls to DMCreateGlobalVector() 1415564755cdSBarry Smith and DMCreateLocalVector() and calls to VecDuplicate() if more are needed. 141647c6ae99SBarry Smith 1417897f7067SBarry Smith You must call DMSetUp() after this call before using this DM. 1418897f7067SBarry Smith 1419897f7067SBarry Smith If you wish to use the options database to change values in the DMDA call DMSetFromOptions() after this call 1420897f7067SBarry Smith but before DMSetUp(). 1421897f7067SBarry Smith 1422db781477SPatrick Sanan .seealso: `DMDestroy()`, `DMView()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMGlobalToLocalBegin()`, `DMDAGetRefinementFactor()`, 1423db781477SPatrick Sanan `DMGlobalToLocalEnd()`, `DMLocalToGlobalBegin()`, `DMLocalToLocalBegin()`, `DMLocalToLocalEnd()`, `DMDASetRefinementFactor()`, 1424db781477SPatrick Sanan `DMDAGetInfo()`, `DMCreateGlobalVector()`, `DMCreateLocalVector()`, `DMDACreateNaturalVector()`, `DMLoad()`, `DMDAGetOwnershipRanges()`, 1425db781477SPatrick Sanan `DMStagCreate3d()` 142647c6ae99SBarry Smith 142747c6ae99SBarry Smith @*/ 1428bff4a2f0SMatthew G. Knepley PetscErrorCode DMDACreate3d(MPI_Comm comm,DMBoundaryType bx,DMBoundaryType by,DMBoundaryType bz,DMDAStencilType stencil_type,PetscInt M, 14299a42bb27SBarry 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) 143047c6ae99SBarry Smith { 143147c6ae99SBarry Smith PetscFunctionBegin; 14329566063dSJacob Faibussowitsch PetscCall(DMDACreate(comm, da)); 14339566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*da, 3)); 14349566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(*da, M, N, P)); 14359566063dSJacob Faibussowitsch PetscCall(DMDASetNumProcs(*da, m, n, p)); 14369566063dSJacob Faibussowitsch PetscCall(DMDASetBoundaryType(*da, bx, by, bz)); 14379566063dSJacob Faibussowitsch PetscCall(DMDASetDof(*da, dof)); 14389566063dSJacob Faibussowitsch PetscCall(DMDASetStencilType(*da, stencil_type)); 14399566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(*da, s)); 14409566063dSJacob Faibussowitsch PetscCall(DMDASetOwnershipRanges(*da, lx, ly, lz)); 144147c6ae99SBarry Smith PetscFunctionReturn(0); 144247c6ae99SBarry Smith } 1443