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)); 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; 8547c6ae99SBarry Smith 869566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(viewer,0,&draw)); 879566063dSJacob Faibussowitsch PetscCall(PetscDrawIsNull(draw,&isnull)); 885b399a63SLisandro Dalcin if (isnull) PetscFunctionReturn(0); 8947c6ae99SBarry Smith 909566063dSJacob Faibussowitsch PetscCall(PetscDrawCheckResizedWindow(draw)); 919566063dSJacob Faibussowitsch PetscCall(PetscDrawClear(draw)); 929566063dSJacob Faibussowitsch PetscCall(PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax)); 935b399a63SLisandro Dalcin 94d0609cedSBarry Smith PetscDrawCollectiveBegin(draw); 9547c6ae99SBarry Smith /* first processor draw all node lines */ 96dd400576SPatrick Sanan if (rank == 0) { 9747c6ae99SBarry Smith for (k=0; k<dd->P; k++) { 9847c6ae99SBarry Smith ymin = 0.0; ymax = (PetscReal)(dd->N - 1); 9947c6ae99SBarry Smith for (xmin=(PetscReal)(k*(dd->M+1)); xmin<(PetscReal)(dd->M+(k*(dd->M+1))); xmin++) { 1009566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK)); 10147c6ae99SBarry Smith } 10247c6ae99SBarry Smith xmin = (PetscReal)(k*(dd->M+1)); xmax = xmin + (PetscReal)(dd->M - 1); 10347c6ae99SBarry Smith for (ymin=0; ymin<(PetscReal)dd->N; ymin++) { 1049566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK)); 10547c6ae99SBarry Smith } 10647c6ae99SBarry Smith } 10747c6ae99SBarry Smith } 108d0609cedSBarry Smith PetscDrawCollectiveEnd(draw); 1099566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 1109566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw)); 11147c6ae99SBarry Smith 112d0609cedSBarry Smith PetscDrawCollectiveBegin(draw); 1135b399a63SLisandro Dalcin /*Go through and draw for each plane*/ 1145b399a63SLisandro Dalcin for (k=0; k<dd->P; k++) { 11547c6ae99SBarry Smith if ((k >= dd->zs) && (k < dd->ze)) { 11647c6ae99SBarry Smith /* draw my box */ 11747c6ae99SBarry Smith ymin = dd->ys; 11847c6ae99SBarry Smith ymax = dd->ye - 1; 11947c6ae99SBarry Smith xmin = dd->xs/dd->w + (dd->M+1)*k; 12047c6ae99SBarry Smith xmax =(dd->xe-1)/dd->w + (dd->M+1)*k; 12147c6ae99SBarry Smith 1229566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED)); 1239566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED)); 1249566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED)); 1259566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED)); 12647c6ae99SBarry Smith 12747c6ae99SBarry Smith xmin = dd->xs/dd->w; 12847c6ae99SBarry Smith xmax =(dd->xe-1)/dd->w; 12947c6ae99SBarry Smith 130832b7cebSLisandro Dalcin /* identify which processor owns the box */ 1319566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(node,sizeof(node),"%d",(int)rank)); 1329566063dSJacob Faibussowitsch PetscCall(PetscDrawString(draw,xmin+(dd->M+1)*k+.2,ymin+.3,PETSC_DRAW_RED,node)); 13347c6ae99SBarry Smith /* put in numbers*/ 13447c6ae99SBarry Smith base = (dd->base+(dd->xe-dd->xs)*(dd->ye-dd->ys)*(k-dd->zs))/dd->w; 13547c6ae99SBarry Smith for (y=ymin; y<=ymax; y++) { 13647c6ae99SBarry Smith for (x=xmin+(dd->M+1)*k; x<=xmax+(dd->M+1)*k; x++) { 1379566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(node,sizeof(node),"%d",(int)base++)); 1389566063dSJacob Faibussowitsch PetscCall(PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node)); 13947c6ae99SBarry Smith } 14047c6ae99SBarry Smith } 14147c6ae99SBarry Smith 14247c6ae99SBarry Smith } 14347c6ae99SBarry Smith } 144d0609cedSBarry Smith PetscDrawCollectiveEnd(draw); 1459566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 1469566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw)); 14747c6ae99SBarry Smith 148d0609cedSBarry Smith PetscDrawCollectiveBegin(draw); 14947c6ae99SBarry Smith for (k=0-dd->s; k<dd->P+dd->s; k++) { 15047c6ae99SBarry Smith /* Go through and draw for each plane */ 15147c6ae99SBarry Smith if ((k >= dd->Zs) && (k < dd->Ze)) { 15247c6ae99SBarry Smith /* overlay ghost numbers, useful for error checking */ 153565245c5SBarry Smith base = (dd->Xe-dd->Xs)*(dd->Ye-dd->Ys)*(k-dd->Zs)/dd->w; 1549566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(da->ltogmap,&idx)); 15547c6ae99SBarry Smith plane=k; 156565245c5SBarry Smith /* Keep z wrap around points on the drawing */ 1578865f1eaSKarl Rupp if (k<0) plane=dd->P+k; 1588865f1eaSKarl Rupp if (k>=dd->P) plane=k-dd->P; 15947c6ae99SBarry Smith ymin = dd->Ys; ymax = dd->Ye; 16047c6ae99SBarry Smith xmin = (dd->M+1)*plane*dd->w; 16147c6ae99SBarry Smith xmax = (dd->M+1)*plane*dd->w+dd->M*dd->w; 16247c6ae99SBarry Smith for (y=ymin; y<ymax; y++) { 16347c6ae99SBarry Smith for (x=xmin+dd->Xs; x<xmin+dd->Xe; x+=dd->w) { 1648ea3bf28SBarry Smith sprintf(node,"%d",(int)(idx[base])); 16547c6ae99SBarry Smith ycoord = y; 16647c6ae99SBarry Smith /*Keep y wrap around points on drawing */ 1678865f1eaSKarl Rupp if (y<0) ycoord = dd->N+y; 1688865f1eaSKarl Rupp if (y>=dd->N) ycoord = y-dd->N; 16947c6ae99SBarry Smith xcoord = x; /* Keep x wrap points on drawing */ 1708865f1eaSKarl Rupp if (x<xmin) xcoord = xmax - (xmin-x); 1718865f1eaSKarl Rupp if (x>=xmax) xcoord = xmin + (x-xmax); 1729566063dSJacob Faibussowitsch PetscCall(PetscDrawString(draw,xcoord/dd->w,ycoord,PETSC_DRAW_BLUE,node)); 173565245c5SBarry Smith base++; 17447c6ae99SBarry Smith } 17547c6ae99SBarry Smith } 1769566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(da->ltogmap,&idx)); 17747c6ae99SBarry Smith } 17847c6ae99SBarry Smith } 179d0609cedSBarry Smith PetscDrawCollectiveEnd(draw); 1809566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 1819566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw)); 1829566063dSJacob Faibussowitsch PetscCall(PetscDrawSave(draw)); 183c9493c35SStefano Zampini } else if (isglvis) { 1849566063dSJacob Faibussowitsch PetscCall(DMView_DA_GLVis(da,viewer)); 1859a42bb27SBarry Smith } else if (isbinary) { 1869566063dSJacob Faibussowitsch PetscCall(DMView_DA_Binary(da,viewer)); 1879a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 1889a42bb27SBarry Smith } else if (ismatlab) { 1899566063dSJacob Faibussowitsch PetscCall(DMView_DA_Matlab(da,viewer)); 1909a42bb27SBarry Smith #endif 19111aeaf0aSBarry Smith } 19247c6ae99SBarry Smith PetscFunctionReturn(0); 19347c6ae99SBarry Smith } 19447c6ae99SBarry Smith 1957087cfbeSBarry Smith PetscErrorCode DMSetUp_DA_3D(DM da) 19647c6ae99SBarry Smith { 19747c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 19847c6ae99SBarry Smith const PetscInt M = dd->M; 19947c6ae99SBarry Smith const PetscInt N = dd->N; 20047c6ae99SBarry Smith const PetscInt P = dd->P; 20147c6ae99SBarry Smith PetscInt m = dd->m; 20247c6ae99SBarry Smith PetscInt n = dd->n; 20347c6ae99SBarry Smith PetscInt p = dd->p; 20447c6ae99SBarry Smith const PetscInt dof = dd->w; 20547c6ae99SBarry Smith const PetscInt s = dd->s; 206bff4a2f0SMatthew G. Knepley DMBoundaryType bx = dd->bx; 207bff4a2f0SMatthew G. Knepley DMBoundaryType by = dd->by; 208bff4a2f0SMatthew G. Knepley DMBoundaryType bz = dd->bz; 20919fd82e9SBarry Smith DMDAStencilType stencil_type = dd->stencil_type; 21047c6ae99SBarry Smith PetscInt *lx = dd->lx; 21147c6ae99SBarry Smith PetscInt *ly = dd->ly; 21247c6ae99SBarry Smith PetscInt *lz = dd->lz; 21347c6ae99SBarry Smith MPI_Comm comm; 21447c6ae99SBarry Smith PetscMPIInt rank,size; 215ce00eea3SSatish Balay PetscInt xs = 0,xe,ys = 0,ye,zs = 0,ze,x = 0,y = 0,z = 0; 216bd1fc5aeSBarry Smith PetscInt Xs,Xe,Ys,Ye,Zs,Ze,IXs,IXe,IYs,IYe,IZs,IZe,pm; 2178ea3bf28SBarry Smith PetscInt left,right,up,down,bottom,top,i,j,k,*idx,nn; 21847c6ae99SBarry Smith PetscInt n0,n1,n2,n3,n4,n5,n6,n7,n8,n9,n10,n11,n12,n14; 21947c6ae99SBarry Smith PetscInt n15,n16,n17,n18,n19,n20,n21,n22,n23,n24,n25,n26; 220db87c5ecSEthan Coon PetscInt *bases,*ldims,base,x_t,y_t,z_t,s_t,count,s_x,s_y,s_z; 22147c6ae99SBarry Smith PetscInt sn0 = 0,sn1 = 0,sn2 = 0,sn3 = 0,sn5 = 0,sn6 = 0,sn7 = 0; 22247c6ae99SBarry Smith PetscInt sn8 = 0,sn9 = 0,sn11 = 0,sn15 = 0,sn24 = 0,sn25 = 0,sn26 = 0; 22347c6ae99SBarry Smith PetscInt sn17 = 0,sn18 = 0,sn19 = 0,sn20 = 0,sn21 = 0,sn23 = 0; 22447c6ae99SBarry Smith Vec local,global; 225bd1fc5aeSBarry Smith VecScatter gtol; 22645b6f7e9SBarry Smith IS to,from; 2276f951b95Secoon PetscBool twod; 22847c6ae99SBarry Smith 22947c6ae99SBarry Smith PetscFunctionBegin; 2307a8be351SBarry 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"); 2319566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject) da, &comm)); 2323855c12bSBarry Smith #if !defined(PETSC_USE_64BIT_INDICES) 23363a3b9bcSJacob 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); 2343855c12bSBarry Smith #endif 2353855c12bSBarry Smith 2369566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 2379566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm,&rank)); 23847c6ae99SBarry Smith 23947c6ae99SBarry Smith if (m != PETSC_DECIDE) { 24063a3b9bcSJacob Faibussowitsch PetscCheck(m >= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %" PetscInt_FMT,m); 24163a3b9bcSJacob Faibussowitsch else PetscCheck(m <= size,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %" PetscInt_FMT " %d",m,size); 24247c6ae99SBarry Smith } 24347c6ae99SBarry Smith if (n != PETSC_DECIDE) { 24463a3b9bcSJacob Faibussowitsch PetscCheck(n >= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %" PetscInt_FMT,n); 24563a3b9bcSJacob Faibussowitsch else PetscCheck(n <= size,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %" PetscInt_FMT " %d",n,size); 24647c6ae99SBarry Smith } 24747c6ae99SBarry Smith if (p != PETSC_DECIDE) { 24863a3b9bcSJacob Faibussowitsch PetscCheck(p >= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Z direction: %" PetscInt_FMT,p); 24963a3b9bcSJacob Faibussowitsch else PetscCheck(p <= size,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Z direction: %" PetscInt_FMT " %d",p,size); 25047c6ae99SBarry Smith } 251*1dca8a05SBarry 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); 25247c6ae99SBarry Smith 25347c6ae99SBarry Smith /* Partition the array among the processors */ 25447c6ae99SBarry Smith if (m == PETSC_DECIDE && n != PETSC_DECIDE && p != PETSC_DECIDE) { 25547c6ae99SBarry Smith m = size/(n*p); 25647c6ae99SBarry Smith } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) { 25747c6ae99SBarry Smith n = size/(m*p); 25847c6ae99SBarry Smith } else if (m != PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) { 25947c6ae99SBarry Smith p = size/(m*n); 26047c6ae99SBarry Smith } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) { 26147c6ae99SBarry Smith /* try for squarish distribution */ 262369cc0aeSBarry Smith m = (int)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)N*p))); 26347c6ae99SBarry Smith if (!m) m = 1; 26447c6ae99SBarry Smith while (m > 0) { 26547c6ae99SBarry Smith n = size/(m*p); 26647c6ae99SBarry Smith if (m*n*p == size) break; 26747c6ae99SBarry Smith m--; 26847c6ae99SBarry Smith } 26963a3b9bcSJacob Faibussowitsch PetscCheck(m,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad p value: p = %" PetscInt_FMT,p); 27047c6ae99SBarry Smith if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;} 27147c6ae99SBarry Smith } else if (m == PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) { 27247c6ae99SBarry Smith /* try for squarish distribution */ 273369cc0aeSBarry Smith m = (int)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)P*n))); 27447c6ae99SBarry Smith if (!m) m = 1; 27547c6ae99SBarry Smith while (m > 0) { 27647c6ae99SBarry Smith p = size/(m*n); 27747c6ae99SBarry Smith if (m*n*p == size) break; 27847c6ae99SBarry Smith m--; 27947c6ae99SBarry Smith } 28063a3b9bcSJacob Faibussowitsch PetscCheck(m,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad n value: n = %" PetscInt_FMT,n); 28147c6ae99SBarry Smith if (M > P && m < p) {PetscInt _m = m; m = p; p = _m;} 28247c6ae99SBarry Smith } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) { 28347c6ae99SBarry Smith /* try for squarish distribution */ 284369cc0aeSBarry Smith n = (int)(0.5 + PetscSqrtReal(((PetscReal)N)*((PetscReal)size)/((PetscReal)P*m))); 28547c6ae99SBarry Smith if (!n) n = 1; 28647c6ae99SBarry Smith while (n > 0) { 28747c6ae99SBarry Smith p = size/(m*n); 28847c6ae99SBarry Smith if (m*n*p == size) break; 28947c6ae99SBarry Smith n--; 29047c6ae99SBarry Smith } 29163a3b9bcSJacob Faibussowitsch PetscCheck(n,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad m value: m = %" PetscInt_FMT,n); 29247c6ae99SBarry Smith if (N > P && n < p) {PetscInt _n = n; n = p; p = _n;} 29347c6ae99SBarry Smith } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) { 29447c6ae99SBarry Smith /* try for squarish distribution */ 295369cc0aeSBarry Smith n = (PetscInt)(0.5 + PetscPowReal(((PetscReal)N*N)*((PetscReal)size)/((PetscReal)P*M),(PetscReal)(1./3.))); 29647c6ae99SBarry Smith if (!n) n = 1; 29747c6ae99SBarry Smith while (n > 0) { 29847c6ae99SBarry Smith pm = size/n; 29947c6ae99SBarry Smith if (n*pm == size) break; 30047c6ae99SBarry Smith n--; 30147c6ae99SBarry Smith } 30247c6ae99SBarry Smith if (!n) n = 1; 303369cc0aeSBarry Smith m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)P*n))); 30447c6ae99SBarry Smith if (!m) m = 1; 30547c6ae99SBarry Smith while (m > 0) { 30647c6ae99SBarry Smith p = size/(m*n); 30747c6ae99SBarry Smith if (m*n*p == size) break; 30847c6ae99SBarry Smith m--; 30947c6ae99SBarry Smith } 31047c6ae99SBarry Smith if (M > P && m < p) {PetscInt _m = m; m = p; p = _m;} 31108401ef6SPierre Jolivet } else PetscCheck(m*n*p == size,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition"); 31247c6ae99SBarry Smith 31308401ef6SPierre Jolivet PetscCheck(m*n*p == size,PetscObjectComm((PetscObject)da),PETSC_ERR_PLIB,"Could not find good partition"); 31463a3b9bcSJacob Faibussowitsch PetscCheck(M >= m,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %" PetscInt_FMT " %" PetscInt_FMT,M,m); 31563a3b9bcSJacob Faibussowitsch PetscCheck(N >= n,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %" PetscInt_FMT " %" PetscInt_FMT,N,n); 31663a3b9bcSJacob Faibussowitsch PetscCheck(P >= p,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Partition in z direction is too fine! %" PetscInt_FMT " %" PetscInt_FMT,P,p); 31747c6ae99SBarry Smith 31847c6ae99SBarry Smith /* 31947c6ae99SBarry Smith Determine locally owned region 32047c6ae99SBarry Smith [x, y, or z]s is the first local node number, [x, y, z] is the number of local nodes 32147c6ae99SBarry Smith */ 32247c6ae99SBarry Smith 32347c6ae99SBarry Smith if (!lx) { 3249566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m, &dd->lx)); 32547c6ae99SBarry Smith lx = dd->lx; 3268865f1eaSKarl Rupp for (i=0; i<m; i++) lx[i] = M/m + ((M % m) > (i % m)); 32747c6ae99SBarry Smith } 32847c6ae99SBarry Smith x = lx[rank % m]; 32947c6ae99SBarry Smith xs = 0; 3308865f1eaSKarl Rupp for (i=0; i<(rank%m); i++) xs += lx[i]; 331*1dca8a05SBarry 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); 33247c6ae99SBarry Smith 33347c6ae99SBarry Smith if (!ly) { 3349566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &dd->ly)); 33547c6ae99SBarry Smith ly = dd->ly; 3368865f1eaSKarl Rupp for (i=0; i<n; i++) ly[i] = N/n + ((N % n) > (i % n)); 33747c6ae99SBarry Smith } 33847c6ae99SBarry Smith y = ly[(rank % (m*n))/m]; 339*1dca8a05SBarry 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); 34030729d88SBarry Smith 34147c6ae99SBarry Smith ys = 0; 3428865f1eaSKarl Rupp for (i=0; i<(rank % (m*n))/m; i++) ys += ly[i]; 34347c6ae99SBarry Smith 34447c6ae99SBarry Smith if (!lz) { 3459566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(p, &dd->lz)); 34647c6ae99SBarry Smith lz = dd->lz; 3478865f1eaSKarl Rupp for (i=0; i<p; i++) lz[i] = P/p + ((P % p) > (i % p)); 34847c6ae99SBarry Smith } 34947c6ae99SBarry Smith z = lz[rank/(m*n)]; 350bcea557cSEthan Coon 351fdc81ce8SEthan Coon /* note this is different than x- and y-, as we will handle as an important special 352bff4a2f0SMatthew G. Knepley case when p=P=1 and DM_BOUNDARY_PERIODIC and s > z. This is to deal with 2D problems 353fdc81ce8SEthan Coon in a 3D code. Additional code for this case is noted with "2d case" comments */ 3546f951b95Secoon twod = PETSC_FALSE; 3558865f1eaSKarl Rupp if (P == 1) twod = PETSC_TRUE; 356*1dca8a05SBarry 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); 35747c6ae99SBarry Smith zs = 0; 3588865f1eaSKarl Rupp for (i=0; i<(rank/(m*n)); i++) zs += lz[i]; 35947c6ae99SBarry Smith ye = ys + y; 36047c6ae99SBarry Smith xe = xs + x; 36147c6ae99SBarry Smith ze = zs + z; 36247c6ae99SBarry Smith 36388661749SPeter Brune /* determine ghost region (Xs) and region scattered into (IXs) */ 364d9c9ebe5SPeter Brune if (xs-s > 0) { 365d9c9ebe5SPeter Brune Xs = xs - s; IXs = xs - s; 36688661749SPeter Brune } else { 3678865f1eaSKarl Rupp if (bx) Xs = xs - s; 3688865f1eaSKarl Rupp else Xs = 0; 36988661749SPeter Brune IXs = 0; 37088661749SPeter Brune } 371d9c9ebe5SPeter Brune if (xe+s <= M) { 372d9c9ebe5SPeter Brune Xe = xe + s; IXe = xe + s; 37388661749SPeter Brune } else { 37488661749SPeter Brune if (bx) { 375d9c9ebe5SPeter Brune Xs = xs - s; Xe = xe + s; 3768865f1eaSKarl Rupp } else Xe = M; 37788661749SPeter Brune IXe = M; 37888661749SPeter Brune } 37947c6ae99SBarry Smith 380bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC || bx == DM_BOUNDARY_MIRROR) { 381d9c9ebe5SPeter Brune IXs = xs - s; 382d9c9ebe5SPeter Brune IXe = xe + s; 383d9c9ebe5SPeter Brune Xs = xs - s; 384d9c9ebe5SPeter Brune Xe = xe + s; 38588661749SPeter Brune } 38688661749SPeter Brune 387d9c9ebe5SPeter Brune if (ys-s > 0) { 388d9c9ebe5SPeter Brune Ys = ys - s; IYs = ys - s; 38988661749SPeter Brune } else { 3908865f1eaSKarl Rupp if (by) Ys = ys - s; 3918865f1eaSKarl Rupp else Ys = 0; 39288661749SPeter Brune IYs = 0; 39388661749SPeter Brune } 394d9c9ebe5SPeter Brune if (ye+s <= N) { 395d9c9ebe5SPeter Brune Ye = ye + s; IYe = ye + s; 39688661749SPeter Brune } else { 3978865f1eaSKarl Rupp if (by) Ye = ye + s; 3988865f1eaSKarl Rupp else Ye = N; 39988661749SPeter Brune IYe = N; 40088661749SPeter Brune } 40188661749SPeter Brune 402bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_PERIODIC || by == DM_BOUNDARY_MIRROR) { 403d9c9ebe5SPeter Brune IYs = ys - s; 404d9c9ebe5SPeter Brune IYe = ye + s; 405d9c9ebe5SPeter Brune Ys = ys - s; 406d9c9ebe5SPeter Brune Ye = ye + s; 40788661749SPeter Brune } 40888661749SPeter Brune 409d9c9ebe5SPeter Brune if (zs-s > 0) { 410d9c9ebe5SPeter Brune Zs = zs - s; IZs = zs - s; 41188661749SPeter Brune } else { 4128865f1eaSKarl Rupp if (bz) Zs = zs - s; 4138865f1eaSKarl Rupp else Zs = 0; 41488661749SPeter Brune IZs = 0; 41588661749SPeter Brune } 416d9c9ebe5SPeter Brune if (ze+s <= P) { 417d9c9ebe5SPeter Brune Ze = ze + s; IZe = ze + s; 41888661749SPeter Brune } else { 4198865f1eaSKarl Rupp if (bz) Ze = ze + s; 4208865f1eaSKarl Rupp else Ze = P; 42188661749SPeter Brune IZe = P; 42288661749SPeter Brune } 42388661749SPeter Brune 424bff4a2f0SMatthew G. Knepley if (bz == DM_BOUNDARY_PERIODIC || bz == DM_BOUNDARY_MIRROR) { 425d9c9ebe5SPeter Brune IZs = zs - s; 426d9c9ebe5SPeter Brune IZe = ze + s; 427d9c9ebe5SPeter Brune Zs = zs - s; 428d9c9ebe5SPeter Brune Ze = ze + s; 42988661749SPeter Brune } 43047c6ae99SBarry Smith 43147c6ae99SBarry Smith /* Resize all X parameters to reflect w */ 432d9c9ebe5SPeter Brune s_x = s; 433d9c9ebe5SPeter Brune s_y = s; 434d9c9ebe5SPeter Brune s_z = s; 43547c6ae99SBarry Smith 43647c6ae99SBarry Smith /* determine starting point of each processor */ 43747c6ae99SBarry Smith nn = x*y*z; 4389566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(size+1,&bases,size,&ldims)); 4399566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm)); 44047c6ae99SBarry Smith bases[0] = 0; 4418865f1eaSKarl Rupp for (i=1; i<=size; i++) bases[i] = ldims[i-1]; 4428865f1eaSKarl Rupp for (i=1; i<=size; i++) bases[i] += bases[i-1]; 443ce00eea3SSatish Balay base = bases[rank]*dof; 44447c6ae99SBarry Smith 44547c6ae99SBarry Smith /* allocate the base parallel and sequential vectors */ 446ce00eea3SSatish Balay dd->Nlocal = x*y*z*dof; 4479566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,NULL,&global)); 448ce00eea3SSatish Balay dd->nlocal = (Xe-Xs)*(Ye-Ys)*(Ze-Zs)*dof; 4499566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,NULL,&local)); 45047c6ae99SBarry Smith 4514104a7a0SPatrick Sanan /* generate global to local vector scatter and local to global mapping*/ 45247c6ae99SBarry Smith 453ce00eea3SSatish Balay /* global to local must include ghost points within the domain, 454ce00eea3SSatish Balay but not ghost points outside the domain that aren't periodic */ 4559566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((IXe-IXs)*(IYe-IYs)*(IZe-IZs),&idx)); 456d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_BOX) { 457ce00eea3SSatish Balay left = IXs - Xs; right = left + (IXe-IXs); 458ce00eea3SSatish Balay bottom = IYs - Ys; top = bottom + (IYe-IYs); 459ce00eea3SSatish Balay down = IZs - Zs; up = down + (IZe-IZs); 460ce00eea3SSatish Balay count = 0; 461ce00eea3SSatish Balay for (i=down; i<up; i++) { 462ce00eea3SSatish Balay for (j=bottom; j<top; j++) { 463ce00eea3SSatish Balay for (k=left; k<right; k++) { 464ce00eea3SSatish Balay idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 465ce00eea3SSatish Balay } 466ce00eea3SSatish Balay } 467ce00eea3SSatish Balay } 4689566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to)); 46947c6ae99SBarry Smith } else { 47047c6ae99SBarry Smith /* This is way ugly! We need to list the funny cross type region */ 471ce00eea3SSatish Balay left = xs - Xs; right = left + x; 47247c6ae99SBarry Smith bottom = ys - Ys; top = bottom + y; 47347c6ae99SBarry Smith down = zs - Zs; up = down + z; 47447c6ae99SBarry Smith count = 0; 475a5b23f4aSJose E. Roman /* the bottom chunk */ 476ce00eea3SSatish Balay for (i=(IZs-Zs); i<down; i++) { 47747c6ae99SBarry Smith for (j=bottom; j<top; j++) { 478ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 47947c6ae99SBarry Smith } 48047c6ae99SBarry Smith } 48147c6ae99SBarry Smith /* the middle piece */ 48247c6ae99SBarry Smith for (i=down; i<up; i++) { 48347c6ae99SBarry Smith /* front */ 484ce00eea3SSatish Balay for (j=(IYs-Ys); j<bottom; j++) { 485ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 48647c6ae99SBarry Smith } 48747c6ae99SBarry Smith /* middle */ 48847c6ae99SBarry Smith for (j=bottom; j<top; j++) { 489ce00eea3SSatish Balay for (k=IXs-Xs; k<IXe-Xs; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 49047c6ae99SBarry Smith } 49147c6ae99SBarry Smith /* back */ 492ce00eea3SSatish Balay for (j=top; j<top+IYe-ye; j++) { 493ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 49447c6ae99SBarry Smith } 49547c6ae99SBarry Smith } 49647c6ae99SBarry Smith /* the top piece */ 497ce00eea3SSatish Balay for (i=up; i<up+IZe-ze; i++) { 49847c6ae99SBarry Smith for (j=bottom; j<top; j++) { 499ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 50047c6ae99SBarry Smith } 50147c6ae99SBarry Smith } 5029566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to)); 50347c6ae99SBarry Smith } 50447c6ae99SBarry Smith 50547c6ae99SBarry Smith /* determine who lies on each side of use stored in n24 n25 n26 50647c6ae99SBarry Smith n21 n22 n23 50747c6ae99SBarry Smith n18 n19 n20 50847c6ae99SBarry Smith 50947c6ae99SBarry Smith n15 n16 n17 51047c6ae99SBarry Smith n12 n14 51147c6ae99SBarry Smith n9 n10 n11 51247c6ae99SBarry Smith 51347c6ae99SBarry Smith n6 n7 n8 51447c6ae99SBarry Smith n3 n4 n5 51547c6ae99SBarry Smith n0 n1 n2 51647c6ae99SBarry Smith */ 51747c6ae99SBarry Smith 51847c6ae99SBarry Smith /* Solve for X,Y, and Z Periodic Case First, Then Modify Solution */ 51947c6ae99SBarry Smith /* Assume Nodes are Internal to the Cube */ 52047c6ae99SBarry Smith n0 = rank - m*n - m - 1; 52147c6ae99SBarry Smith n1 = rank - m*n - m; 52247c6ae99SBarry Smith n2 = rank - m*n - m + 1; 52347c6ae99SBarry Smith n3 = rank - m*n -1; 52447c6ae99SBarry Smith n4 = rank - m*n; 52547c6ae99SBarry Smith n5 = rank - m*n + 1; 52647c6ae99SBarry Smith n6 = rank - m*n + m - 1; 52747c6ae99SBarry Smith n7 = rank - m*n + m; 52847c6ae99SBarry Smith n8 = rank - m*n + m + 1; 52947c6ae99SBarry Smith 53047c6ae99SBarry Smith n9 = rank - m - 1; 53147c6ae99SBarry Smith n10 = rank - m; 53247c6ae99SBarry Smith n11 = rank - m + 1; 53347c6ae99SBarry Smith n12 = rank - 1; 53447c6ae99SBarry Smith n14 = rank + 1; 53547c6ae99SBarry Smith n15 = rank + m - 1; 53647c6ae99SBarry Smith n16 = rank + m; 53747c6ae99SBarry Smith n17 = rank + m + 1; 53847c6ae99SBarry Smith 53947c6ae99SBarry Smith n18 = rank + m*n - m - 1; 54047c6ae99SBarry Smith n19 = rank + m*n - m; 54147c6ae99SBarry Smith n20 = rank + m*n - m + 1; 54247c6ae99SBarry Smith n21 = rank + m*n - 1; 54347c6ae99SBarry Smith n22 = rank + m*n; 54447c6ae99SBarry Smith n23 = rank + m*n + 1; 54547c6ae99SBarry Smith n24 = rank + m*n + m - 1; 54647c6ae99SBarry Smith n25 = rank + m*n + m; 54747c6ae99SBarry Smith n26 = rank + m*n + m + 1; 54847c6ae99SBarry Smith 54947c6ae99SBarry Smith /* Assume Pieces are on Faces of Cube */ 55047c6ae99SBarry Smith 55147c6ae99SBarry Smith if (xs == 0) { /* First assume not corner or edge */ 55247c6ae99SBarry Smith n0 = rank -1 - (m*n); 55347c6ae99SBarry Smith n3 = rank + m -1 - (m*n); 55447c6ae99SBarry Smith n6 = rank + 2*m -1 - (m*n); 55547c6ae99SBarry Smith n9 = rank -1; 55647c6ae99SBarry Smith n12 = rank + m -1; 55747c6ae99SBarry Smith n15 = rank + 2*m -1; 55847c6ae99SBarry Smith n18 = rank -1 + (m*n); 55947c6ae99SBarry Smith n21 = rank + m -1 + (m*n); 56047c6ae99SBarry Smith n24 = rank + 2*m -1 + (m*n); 56147c6ae99SBarry Smith } 56247c6ae99SBarry Smith 563ce00eea3SSatish Balay if (xe == M) { /* First assume not corner or edge */ 56447c6ae99SBarry Smith n2 = rank -2*m +1 - (m*n); 56547c6ae99SBarry Smith n5 = rank - m +1 - (m*n); 56647c6ae99SBarry Smith n8 = rank +1 - (m*n); 56747c6ae99SBarry Smith n11 = rank -2*m +1; 56847c6ae99SBarry Smith n14 = rank - m +1; 56947c6ae99SBarry Smith n17 = rank +1; 57047c6ae99SBarry Smith n20 = rank -2*m +1 + (m*n); 57147c6ae99SBarry Smith n23 = rank - m +1 + (m*n); 57247c6ae99SBarry Smith n26 = rank +1 + (m*n); 57347c6ae99SBarry Smith } 57447c6ae99SBarry Smith 57547c6ae99SBarry Smith if (ys==0) { /* First assume not corner or edge */ 57647c6ae99SBarry Smith n0 = rank + m * (n-1) -1 - (m*n); 57747c6ae99SBarry Smith n1 = rank + m * (n-1) - (m*n); 57847c6ae99SBarry Smith n2 = rank + m * (n-1) +1 - (m*n); 57947c6ae99SBarry Smith n9 = rank + m * (n-1) -1; 58047c6ae99SBarry Smith n10 = rank + m * (n-1); 58147c6ae99SBarry Smith n11 = rank + m * (n-1) +1; 58247c6ae99SBarry Smith n18 = rank + m * (n-1) -1 + (m*n); 58347c6ae99SBarry Smith n19 = rank + m * (n-1) + (m*n); 58447c6ae99SBarry Smith n20 = rank + m * (n-1) +1 + (m*n); 58547c6ae99SBarry Smith } 58647c6ae99SBarry Smith 58747c6ae99SBarry Smith if (ye == N) { /* First assume not corner or edge */ 58847c6ae99SBarry Smith n6 = rank - m * (n-1) -1 - (m*n); 58947c6ae99SBarry Smith n7 = rank - m * (n-1) - (m*n); 59047c6ae99SBarry Smith n8 = rank - m * (n-1) +1 - (m*n); 59147c6ae99SBarry Smith n15 = rank - m * (n-1) -1; 59247c6ae99SBarry Smith n16 = rank - m * (n-1); 59347c6ae99SBarry Smith n17 = rank - m * (n-1) +1; 59447c6ae99SBarry Smith n24 = rank - m * (n-1) -1 + (m*n); 59547c6ae99SBarry Smith n25 = rank - m * (n-1) + (m*n); 59647c6ae99SBarry Smith n26 = rank - m * (n-1) +1 + (m*n); 59747c6ae99SBarry Smith } 59847c6ae99SBarry Smith 59947c6ae99SBarry Smith if (zs == 0) { /* First assume not corner or edge */ 60047c6ae99SBarry Smith n0 = size - (m*n) + rank - m - 1; 60147c6ae99SBarry Smith n1 = size - (m*n) + rank - m; 60247c6ae99SBarry Smith n2 = size - (m*n) + rank - m + 1; 60347c6ae99SBarry Smith n3 = size - (m*n) + rank - 1; 60447c6ae99SBarry Smith n4 = size - (m*n) + rank; 60547c6ae99SBarry Smith n5 = size - (m*n) + rank + 1; 60647c6ae99SBarry Smith n6 = size - (m*n) + rank + m - 1; 60747c6ae99SBarry Smith n7 = size - (m*n) + rank + m; 60847c6ae99SBarry Smith n8 = size - (m*n) + rank + m + 1; 60947c6ae99SBarry Smith } 61047c6ae99SBarry Smith 61147c6ae99SBarry Smith if (ze == P) { /* First assume not corner or edge */ 61247c6ae99SBarry Smith n18 = (m*n) - (size-rank) - m - 1; 61347c6ae99SBarry Smith n19 = (m*n) - (size-rank) - m; 61447c6ae99SBarry Smith n20 = (m*n) - (size-rank) - m + 1; 61547c6ae99SBarry Smith n21 = (m*n) - (size-rank) - 1; 61647c6ae99SBarry Smith n22 = (m*n) - (size-rank); 61747c6ae99SBarry Smith n23 = (m*n) - (size-rank) + 1; 61847c6ae99SBarry Smith n24 = (m*n) - (size-rank) + m - 1; 61947c6ae99SBarry Smith n25 = (m*n) - (size-rank) + m; 62047c6ae99SBarry Smith n26 = (m*n) - (size-rank) + m + 1; 62147c6ae99SBarry Smith } 62247c6ae99SBarry Smith 62347c6ae99SBarry Smith if ((xs==0) && (zs==0)) { /* Assume an edge, not corner */ 62447c6ae99SBarry Smith n0 = size - m*n + rank + m-1 - m; 62547c6ae99SBarry Smith n3 = size - m*n + rank + m-1; 62647c6ae99SBarry Smith n6 = size - m*n + rank + m-1 + m; 62747c6ae99SBarry Smith } 62847c6ae99SBarry Smith 62947c6ae99SBarry Smith if ((xs==0) && (ze==P)) { /* Assume an edge, not corner */ 63047c6ae99SBarry Smith n18 = m*n - (size - rank) + m-1 - m; 63147c6ae99SBarry Smith n21 = m*n - (size - rank) + m-1; 63247c6ae99SBarry Smith n24 = m*n - (size - rank) + m-1 + m; 63347c6ae99SBarry Smith } 63447c6ae99SBarry Smith 63547c6ae99SBarry Smith if ((xs==0) && (ys==0)) { /* Assume an edge, not corner */ 63647c6ae99SBarry Smith n0 = rank + m*n -1 - m*n; 63747c6ae99SBarry Smith n9 = rank + m*n -1; 63847c6ae99SBarry Smith n18 = rank + m*n -1 + m*n; 63947c6ae99SBarry Smith } 64047c6ae99SBarry Smith 64147c6ae99SBarry Smith if ((xs==0) && (ye==N)) { /* Assume an edge, not corner */ 64247c6ae99SBarry Smith n6 = rank - m*(n-1) + m-1 - m*n; 64347c6ae99SBarry Smith n15 = rank - m*(n-1) + m-1; 64447c6ae99SBarry Smith n24 = rank - m*(n-1) + m-1 + m*n; 64547c6ae99SBarry Smith } 64647c6ae99SBarry Smith 647ce00eea3SSatish Balay if ((xe==M) && (zs==0)) { /* Assume an edge, not corner */ 64847c6ae99SBarry Smith n2 = size - (m*n-rank) - (m-1) - m; 64947c6ae99SBarry Smith n5 = size - (m*n-rank) - (m-1); 65047c6ae99SBarry Smith n8 = size - (m*n-rank) - (m-1) + m; 65147c6ae99SBarry Smith } 65247c6ae99SBarry Smith 653ce00eea3SSatish Balay if ((xe==M) && (ze==P)) { /* Assume an edge, not corner */ 65447c6ae99SBarry Smith n20 = m*n - (size - rank) - (m-1) - m; 65547c6ae99SBarry Smith n23 = m*n - (size - rank) - (m-1); 65647c6ae99SBarry Smith n26 = m*n - (size - rank) - (m-1) + m; 65747c6ae99SBarry Smith } 65847c6ae99SBarry Smith 659ce00eea3SSatish Balay if ((xe==M) && (ys==0)) { /* Assume an edge, not corner */ 66047c6ae99SBarry Smith n2 = rank + m*(n-1) - (m-1) - m*n; 66147c6ae99SBarry Smith n11 = rank + m*(n-1) - (m-1); 66247c6ae99SBarry Smith n20 = rank + m*(n-1) - (m-1) + m*n; 66347c6ae99SBarry Smith } 66447c6ae99SBarry Smith 665ce00eea3SSatish Balay if ((xe==M) && (ye==N)) { /* Assume an edge, not corner */ 66647c6ae99SBarry Smith n8 = rank - m*n +1 - m*n; 66747c6ae99SBarry Smith n17 = rank - m*n +1; 66847c6ae99SBarry Smith n26 = rank - m*n +1 + m*n; 66947c6ae99SBarry Smith } 67047c6ae99SBarry Smith 67147c6ae99SBarry Smith if ((ys==0) && (zs==0)) { /* Assume an edge, not corner */ 67247c6ae99SBarry Smith n0 = size - m + rank -1; 67347c6ae99SBarry Smith n1 = size - m + rank; 67447c6ae99SBarry Smith n2 = size - m + rank +1; 67547c6ae99SBarry Smith } 67647c6ae99SBarry Smith 67747c6ae99SBarry Smith if ((ys==0) && (ze==P)) { /* Assume an edge, not corner */ 67847c6ae99SBarry Smith n18 = m*n - (size - rank) + m*(n-1) -1; 67947c6ae99SBarry Smith n19 = m*n - (size - rank) + m*(n-1); 68047c6ae99SBarry Smith n20 = m*n - (size - rank) + m*(n-1) +1; 68147c6ae99SBarry Smith } 68247c6ae99SBarry Smith 68347c6ae99SBarry Smith if ((ye==N) && (zs==0)) { /* Assume an edge, not corner */ 68447c6ae99SBarry Smith n6 = size - (m*n-rank) - m * (n-1) -1; 68547c6ae99SBarry Smith n7 = size - (m*n-rank) - m * (n-1); 68647c6ae99SBarry Smith n8 = size - (m*n-rank) - m * (n-1) +1; 68747c6ae99SBarry Smith } 68847c6ae99SBarry Smith 68947c6ae99SBarry Smith if ((ye==N) && (ze==P)) { /* Assume an edge, not corner */ 69047c6ae99SBarry Smith n24 = rank - (size-m) -1; 69147c6ae99SBarry Smith n25 = rank - (size-m); 69247c6ae99SBarry Smith n26 = rank - (size-m) +1; 69347c6ae99SBarry Smith } 69447c6ae99SBarry Smith 69547c6ae99SBarry Smith /* Check for Corners */ 6968865f1eaSKarl Rupp if ((xs==0) && (ys==0) && (zs==0)) n0 = size -1; 6978865f1eaSKarl Rupp if ((xs==0) && (ys==0) && (ze==P)) n18 = m*n-1; 6988865f1eaSKarl Rupp if ((xs==0) && (ye==N) && (zs==0)) n6 = (size-1)-m*(n-1); 6998865f1eaSKarl Rupp if ((xs==0) && (ye==N) && (ze==P)) n24 = m-1; 7008865f1eaSKarl Rupp if ((xe==M) && (ys==0) && (zs==0)) n2 = size-m; 7018865f1eaSKarl Rupp if ((xe==M) && (ys==0) && (ze==P)) n20 = m*n-m; 7028865f1eaSKarl Rupp if ((xe==M) && (ye==N) && (zs==0)) n8 = size-m*n; 7038865f1eaSKarl Rupp if ((xe==M) && (ye==N) && (ze==P)) n26 = 0; 70447c6ae99SBarry Smith 70547c6ae99SBarry Smith /* Check for when not X,Y, and Z Periodic */ 70647c6ae99SBarry Smith 70747c6ae99SBarry Smith /* If not X periodic */ 708bff4a2f0SMatthew G. Knepley if (bx != DM_BOUNDARY_PERIODIC) { 7098865f1eaSKarl Rupp if (xs==0) n0 = n3 = n6 = n9 = n12 = n15 = n18 = n21 = n24 = -2; 7108865f1eaSKarl Rupp if (xe==M) n2 = n5 = n8 = n11 = n14 = n17 = n20 = n23 = n26 = -2; 71147c6ae99SBarry Smith } 71247c6ae99SBarry Smith 71347c6ae99SBarry Smith /* If not Y periodic */ 714bff4a2f0SMatthew G. Knepley if (by != DM_BOUNDARY_PERIODIC) { 7158865f1eaSKarl Rupp if (ys==0) n0 = n1 = n2 = n9 = n10 = n11 = n18 = n19 = n20 = -2; 7168865f1eaSKarl Rupp if (ye==N) n6 = n7 = n8 = n15 = n16 = n17 = n24 = n25 = n26 = -2; 71747c6ae99SBarry Smith } 71847c6ae99SBarry Smith 71947c6ae99SBarry Smith /* If not Z periodic */ 720bff4a2f0SMatthew G. Knepley if (bz != DM_BOUNDARY_PERIODIC) { 7218865f1eaSKarl Rupp if (zs==0) n0 = n1 = n2 = n3 = n4 = n5 = n6 = n7 = n8 = -2; 7228865f1eaSKarl Rupp if (ze==P) n18 = n19 = n20 = n21 = n22 = n23 = n24 = n25 = n26 = -2; 72347c6ae99SBarry Smith } 72447c6ae99SBarry Smith 7259566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(27,&dd->neighbors)); 7268865f1eaSKarl Rupp 72747c6ae99SBarry Smith dd->neighbors[0] = n0; 72847c6ae99SBarry Smith dd->neighbors[1] = n1; 72947c6ae99SBarry Smith dd->neighbors[2] = n2; 73047c6ae99SBarry Smith dd->neighbors[3] = n3; 73147c6ae99SBarry Smith dd->neighbors[4] = n4; 73247c6ae99SBarry Smith dd->neighbors[5] = n5; 73347c6ae99SBarry Smith dd->neighbors[6] = n6; 73447c6ae99SBarry Smith dd->neighbors[7] = n7; 73547c6ae99SBarry Smith dd->neighbors[8] = n8; 73647c6ae99SBarry Smith dd->neighbors[9] = n9; 73747c6ae99SBarry Smith dd->neighbors[10] = n10; 73847c6ae99SBarry Smith dd->neighbors[11] = n11; 73947c6ae99SBarry Smith dd->neighbors[12] = n12; 74047c6ae99SBarry Smith dd->neighbors[13] = rank; 74147c6ae99SBarry Smith dd->neighbors[14] = n14; 74247c6ae99SBarry Smith dd->neighbors[15] = n15; 74347c6ae99SBarry Smith dd->neighbors[16] = n16; 74447c6ae99SBarry Smith dd->neighbors[17] = n17; 74547c6ae99SBarry Smith dd->neighbors[18] = n18; 74647c6ae99SBarry Smith dd->neighbors[19] = n19; 74747c6ae99SBarry Smith dd->neighbors[20] = n20; 74847c6ae99SBarry Smith dd->neighbors[21] = n21; 74947c6ae99SBarry Smith dd->neighbors[22] = n22; 75047c6ae99SBarry Smith dd->neighbors[23] = n23; 75147c6ae99SBarry Smith dd->neighbors[24] = n24; 75247c6ae99SBarry Smith dd->neighbors[25] = n25; 75347c6ae99SBarry Smith dd->neighbors[26] = n26; 75447c6ae99SBarry Smith 75547c6ae99SBarry Smith /* If star stencil then delete the corner neighbors */ 756d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_STAR) { 75747c6ae99SBarry Smith /* save information about corner neighbors */ 75847c6ae99SBarry Smith sn0 = n0; sn1 = n1; sn2 = n2; sn3 = n3; sn5 = n5; sn6 = n6; sn7 = n7; 75947c6ae99SBarry Smith sn8 = n8; sn9 = n9; sn11 = n11; sn15 = n15; sn17 = n17; sn18 = n18; 76047c6ae99SBarry Smith sn19 = n19; sn20 = n20; sn21 = n21; sn23 = n23; sn24 = n24; sn25 = n25; 76147c6ae99SBarry Smith sn26 = n26; 7628865f1eaSKarl Rupp n0 = n1 = n2 = n3 = n5 = n6 = n7 = n8 = n9 = n11 = n15 = n17 = n18 = n19 = n20 = n21 = n23 = n24 = n25 = n26 = -1; 76347c6ae99SBarry Smith } 76447c6ae99SBarry Smith 7659566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((Xe-Xs)*(Ye-Ys)*(Ze-Zs),&idx)); 76647c6ae99SBarry Smith 76747c6ae99SBarry Smith nn = 0; 76847c6ae99SBarry Smith /* Bottom Level */ 76947c6ae99SBarry Smith for (k=0; k<s_z; k++) { 77047c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 77147c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 772ce00eea3SSatish Balay x_t = lx[n0 % m]; 77347c6ae99SBarry Smith y_t = ly[(n0 % (m*n))/m]; 77447c6ae99SBarry Smith z_t = lz[n0 / (m*n)]; 77547c6ae99SBarry 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; 7768865f1eaSKarl 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 */ 7778865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 77847c6ae99SBarry Smith } 77947c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 78047c6ae99SBarry Smith x_t = x; 78147c6ae99SBarry Smith y_t = ly[(n1 % (m*n))/m]; 78247c6ae99SBarry Smith z_t = lz[n1 / (m*n)]; 78347c6ae99SBarry 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; 7848865f1eaSKarl Rupp if (twod && (s_t < 0)) s_t = bases[n1] + x_t*y_t*z_t - (s_y+1-i)*x_t; /* 2D case */ 7858865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 78647c6ae99SBarry Smith } 78747c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 788ce00eea3SSatish Balay x_t = lx[n2 % m]; 78947c6ae99SBarry Smith y_t = ly[(n2 % (m*n))/m]; 79047c6ae99SBarry Smith z_t = lz[n2 / (m*n)]; 79147c6ae99SBarry 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; 7928865f1eaSKarl Rupp if (twod && (s_t < 0)) s_t = bases[n2] + x_t*y_t*z_t - (s_y+1-i)*x_t; /* 2D case */ 7938865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 79447c6ae99SBarry Smith } 79547c6ae99SBarry Smith } 79647c6ae99SBarry Smith 79747c6ae99SBarry Smith for (i=0; i<y; i++) { 79847c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 799ce00eea3SSatish Balay x_t = lx[n3 % m]; 80047c6ae99SBarry Smith y_t = y; 80147c6ae99SBarry Smith z_t = lz[n3 / (m*n)]; 80247c6ae99SBarry Smith s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8038865f1eaSKarl 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 */ 8048865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 80547c6ae99SBarry Smith } 80647c6ae99SBarry Smith 80747c6ae99SBarry Smith if (n4 >= 0) { /* middle */ 80847c6ae99SBarry Smith x_t = x; 80947c6ae99SBarry Smith y_t = y; 81047c6ae99SBarry Smith z_t = lz[n4 / (m*n)]; 81147c6ae99SBarry Smith s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8128865f1eaSKarl 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 */ 8138865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 814bff4a2f0SMatthew G. Knepley } else if (bz == DM_BOUNDARY_MIRROR) { 815a6fd6b0aSBarry Smith for (j=0; j<x; j++) idx[nn++] = bases[rank] + j + x*i + (s_z - k - 1)*x*y; 81647c6ae99SBarry Smith } 81747c6ae99SBarry Smith 81847c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 819ce00eea3SSatish Balay x_t = lx[n5 % m]; 82047c6ae99SBarry Smith y_t = y; 82147c6ae99SBarry Smith z_t = lz[n5 / (m*n)]; 82247c6ae99SBarry Smith s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8238865f1eaSKarl 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 */ 8248865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 82547c6ae99SBarry Smith } 82647c6ae99SBarry Smith } 82747c6ae99SBarry Smith 82847c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 82947c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 830ce00eea3SSatish Balay x_t = lx[n6 % m]; 83147c6ae99SBarry Smith y_t = ly[(n6 % (m*n))/m]; 83247c6ae99SBarry Smith z_t = lz[n6 / (m*n)]; 83347c6ae99SBarry Smith s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8348865f1eaSKarl 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 */ 8358865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 83647c6ae99SBarry Smith } 83747c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 83847c6ae99SBarry Smith x_t = x; 83947c6ae99SBarry Smith y_t = ly[(n7 % (m*n))/m]; 84047c6ae99SBarry Smith z_t = lz[n7 / (m*n)]; 84147c6ae99SBarry Smith s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8428865f1eaSKarl 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 */ 8438865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 84447c6ae99SBarry Smith } 84547c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 846ce00eea3SSatish Balay x_t = lx[n8 % m]; 84747c6ae99SBarry Smith y_t = ly[(n8 % (m*n))/m]; 84847c6ae99SBarry Smith z_t = lz[n8 / (m*n)]; 84947c6ae99SBarry Smith s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8508865f1eaSKarl 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 */ 8518865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 85247c6ae99SBarry Smith } 85347c6ae99SBarry Smith } 85447c6ae99SBarry Smith } 85547c6ae99SBarry Smith 85647c6ae99SBarry Smith /* Middle Level */ 85747c6ae99SBarry Smith for (k=0; k<z; k++) { 85847c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 85947c6ae99SBarry Smith if (n9 >= 0) { /* left below */ 860ce00eea3SSatish Balay x_t = lx[n9 % m]; 86147c6ae99SBarry Smith y_t = ly[(n9 % (m*n))/m]; 86247c6ae99SBarry Smith /* z_t = z; */ 86347c6ae99SBarry Smith s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 8648865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 86547c6ae99SBarry Smith } 86647c6ae99SBarry Smith if (n10 >= 0) { /* directly below */ 86747c6ae99SBarry Smith x_t = x; 86847c6ae99SBarry Smith y_t = ly[(n10 % (m*n))/m]; 86947c6ae99SBarry Smith /* z_t = z; */ 87047c6ae99SBarry Smith s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 8718865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 872bff4a2f0SMatthew G. Knepley } else if (by == DM_BOUNDARY_MIRROR) { 873a6fd6b0aSBarry Smith for (j=0; j<x; j++) idx[nn++] = bases[rank] + j + k*x*y + (s_y - i)*x; 87447c6ae99SBarry Smith } 87547c6ae99SBarry Smith if (n11 >= 0) { /* right below */ 876ce00eea3SSatish Balay x_t = lx[n11 % m]; 87747c6ae99SBarry Smith y_t = ly[(n11 % (m*n))/m]; 87847c6ae99SBarry Smith /* z_t = z; */ 87947c6ae99SBarry Smith s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 8808865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 88147c6ae99SBarry Smith } 88247c6ae99SBarry Smith } 88347c6ae99SBarry Smith 88447c6ae99SBarry Smith for (i=0; i<y; i++) { 88547c6ae99SBarry Smith if (n12 >= 0) { /* directly left */ 886ce00eea3SSatish Balay x_t = lx[n12 % m]; 88747c6ae99SBarry Smith y_t = y; 88847c6ae99SBarry Smith /* z_t = z; */ 88947c6ae99SBarry Smith s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t; 8908865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 891bff4a2f0SMatthew G. Knepley } else if (bx == DM_BOUNDARY_MIRROR) { 892a6fd6b0aSBarry Smith for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + s_x - j - 1 + k*x*y + i*x; 89347c6ae99SBarry Smith } 89447c6ae99SBarry Smith 89547c6ae99SBarry Smith /* Interior */ 89647c6ae99SBarry Smith s_t = bases[rank] + i*x + k*x*y; 8978865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = s_t++; 89847c6ae99SBarry Smith 89947c6ae99SBarry Smith if (n14 >= 0) { /* directly right */ 900ce00eea3SSatish Balay x_t = lx[n14 % m]; 90147c6ae99SBarry Smith y_t = y; 90247c6ae99SBarry Smith /* z_t = z; */ 90347c6ae99SBarry Smith s_t = bases[n14] + i*x_t + k*x_t*y_t; 9048865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 905bff4a2f0SMatthew G. Knepley } else if (bx == DM_BOUNDARY_MIRROR) { 906a6fd6b0aSBarry Smith for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x - j - 1 + k*x*y + i*x; 90747c6ae99SBarry Smith } 90847c6ae99SBarry Smith } 90947c6ae99SBarry Smith 91047c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 91147c6ae99SBarry Smith if (n15 >= 0) { /* left above */ 912ce00eea3SSatish Balay x_t = lx[n15 % m]; 91347c6ae99SBarry Smith y_t = ly[(n15 % (m*n))/m]; 91447c6ae99SBarry Smith /* z_t = z; */ 91547c6ae99SBarry Smith s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t; 9168865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 91747c6ae99SBarry Smith } 91847c6ae99SBarry Smith if (n16 >= 0) { /* directly above */ 91947c6ae99SBarry Smith x_t = x; 92047c6ae99SBarry Smith y_t = ly[(n16 % (m*n))/m]; 92147c6ae99SBarry Smith /* z_t = z; */ 92247c6ae99SBarry Smith s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t; 9238865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 924bff4a2f0SMatthew G. Knepley } else if (by == DM_BOUNDARY_MIRROR) { 925a6fd6b0aSBarry Smith for (j=0; j<x; j++) idx[nn++] = bases[rank] + j + k*x*y + (y-i)*x; 92647c6ae99SBarry Smith } 92747c6ae99SBarry Smith if (n17 >= 0) { /* right above */ 928ce00eea3SSatish Balay x_t = lx[n17 % m]; 92947c6ae99SBarry Smith y_t = ly[(n17 % (m*n))/m]; 93047c6ae99SBarry Smith /* z_t = z; */ 93147c6ae99SBarry Smith s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t; 9328865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 93347c6ae99SBarry Smith } 93447c6ae99SBarry Smith } 93547c6ae99SBarry Smith } 93647c6ae99SBarry Smith 93747c6ae99SBarry Smith /* Upper Level */ 93847c6ae99SBarry Smith for (k=0; k<s_z; k++) { 93947c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 94047c6ae99SBarry Smith if (n18 >= 0) { /* left below */ 941ce00eea3SSatish Balay x_t = lx[n18 % m]; 94247c6ae99SBarry Smith y_t = ly[(n18 % (m*n))/m]; 94347c6ae99SBarry Smith /* z_t = lz[n18 / (m*n)]; */ 94447c6ae99SBarry Smith s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 9458865f1eaSKarl 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 */ 9468865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 94747c6ae99SBarry Smith } 94847c6ae99SBarry Smith if (n19 >= 0) { /* directly below */ 94947c6ae99SBarry Smith x_t = x; 95047c6ae99SBarry Smith y_t = ly[(n19 % (m*n))/m]; 95147c6ae99SBarry Smith /* z_t = lz[n19 / (m*n)]; */ 95247c6ae99SBarry Smith s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 9538865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n19] - (s_y+1-i)*x_t + x_t*y_t; /* 2d case */ 9548865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 95547c6ae99SBarry Smith } 95647c6ae99SBarry Smith if (n20 >= 0) { /* right below */ 957ce00eea3SSatish Balay x_t = lx[n20 % m]; 95847c6ae99SBarry Smith y_t = ly[(n20 % (m*n))/m]; 95947c6ae99SBarry Smith /* z_t = lz[n20 / (m*n)]; */ 96047c6ae99SBarry Smith s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 9618865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n20] - (s_y+1-i)*x_t + x_t*y_t; /* 2d case */ 9628865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 96347c6ae99SBarry Smith } 96447c6ae99SBarry Smith } 96547c6ae99SBarry Smith 96647c6ae99SBarry Smith for (i=0; i<y; i++) { 96747c6ae99SBarry Smith if (n21 >= 0) { /* directly left */ 968ce00eea3SSatish Balay x_t = lx[n21 % m]; 96947c6ae99SBarry Smith y_t = y; 97047c6ae99SBarry Smith /* z_t = lz[n21 / (m*n)]; */ 97147c6ae99SBarry Smith s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t; 9728865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n21] + (i+1)*x_t - s_x; /* 2d case */ 9738865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 97447c6ae99SBarry Smith } 97547c6ae99SBarry Smith 97647c6ae99SBarry Smith if (n22 >= 0) { /* middle */ 97747c6ae99SBarry Smith x_t = x; 97847c6ae99SBarry Smith y_t = y; 97947c6ae99SBarry Smith /* z_t = lz[n22 / (m*n)]; */ 98047c6ae99SBarry Smith s_t = bases[n22] + i*x_t + k*x_t*y_t; 9818865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n22] + i*x_t; /* 2d case */ 9828865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 983bff4a2f0SMatthew G. Knepley } else if (bz == DM_BOUNDARY_MIRROR) { 984a6fd6b0aSBarry Smith for (j=0; j<x; j++) idx[nn++] = bases[rank] + j + (z-k-1)*x*y + i*x; 98547c6ae99SBarry Smith } 98647c6ae99SBarry Smith 98747c6ae99SBarry Smith if (n23 >= 0) { /* directly right */ 988ce00eea3SSatish Balay x_t = lx[n23 % m]; 98947c6ae99SBarry Smith y_t = y; 99047c6ae99SBarry Smith /* z_t = lz[n23 / (m*n)]; */ 99147c6ae99SBarry Smith s_t = bases[n23] + i*x_t + k*x_t*y_t; 9928865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n23] + i*x_t; /* 2d case */ 9938865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 99447c6ae99SBarry Smith } 99547c6ae99SBarry Smith } 99647c6ae99SBarry Smith 99747c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 99847c6ae99SBarry Smith if (n24 >= 0) { /* left above */ 999ce00eea3SSatish Balay x_t = lx[n24 % m]; 100047c6ae99SBarry Smith y_t = ly[(n24 % (m*n))/m]; 100147c6ae99SBarry Smith /* z_t = lz[n24 / (m*n)]; */ 100247c6ae99SBarry Smith s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t; 10038865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n24] + i*x_t - s_x; /* 2d case */ 10048865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 100547c6ae99SBarry Smith } 100647c6ae99SBarry Smith if (n25 >= 0) { /* directly above */ 100747c6ae99SBarry Smith x_t = x; 100847c6ae99SBarry Smith y_t = ly[(n25 % (m*n))/m]; 100947c6ae99SBarry Smith /* z_t = lz[n25 / (m*n)]; */ 101047c6ae99SBarry Smith s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t; 10118865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n25] + (i-1)*x_t; /* 2d case */ 10128865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 101347c6ae99SBarry Smith } 101447c6ae99SBarry Smith if (n26 >= 0) { /* right above */ 1015ce00eea3SSatish Balay x_t = lx[n26 % m]; 101647c6ae99SBarry Smith y_t = ly[(n26 % (m*n))/m]; 101747c6ae99SBarry Smith /* z_t = lz[n26 / (m*n)]; */ 101847c6ae99SBarry Smith s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t; 10198865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n26] + (i-1)*x_t; /* 2d case */ 10208865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 102147c6ae99SBarry Smith } 102247c6ae99SBarry Smith } 102347c6ae99SBarry Smith } 1024ce00eea3SSatish Balay 10259566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(comm,dof,nn,idx,PETSC_USE_POINTER,&from)); 10269566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(global,from,local,to,>ol)); 10279566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)da,(PetscObject)gtol)); 10289566063dSJacob Faibussowitsch PetscCall(ISDestroy(&to)); 10299566063dSJacob Faibussowitsch PetscCall(ISDestroy(&from)); 103047c6ae99SBarry Smith 1031d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_STAR) { 103247c6ae99SBarry Smith n0 = sn0; n1 = sn1; n2 = sn2; n3 = sn3; n5 = sn5; n6 = sn6; n7 = sn7; 103347c6ae99SBarry Smith n8 = sn8; n9 = sn9; n11 = sn11; n15 = sn15; n17 = sn17; n18 = sn18; 103447c6ae99SBarry Smith n19 = sn19; n20 = sn20; n21 = sn21; n23 = sn23; n24 = sn24; n25 = sn25; 103547c6ae99SBarry Smith n26 = sn26; 1036ce00eea3SSatish Balay } 103747c6ae99SBarry Smith 1038288e7d53SBarry Smith if (((stencil_type == DMDA_STENCIL_STAR) || (bx != DM_BOUNDARY_PERIODIC && bx) || (by != DM_BOUNDARY_PERIODIC && by) || (bz != DM_BOUNDARY_PERIODIC && bz))) { 1039ce00eea3SSatish Balay /* 1040ce00eea3SSatish Balay Recompute the local to global mappings, this time keeping the 1041ce00eea3SSatish Balay information about the cross corner processor numbers. 1042ce00eea3SSatish Balay */ 104347c6ae99SBarry Smith nn = 0; 104447c6ae99SBarry Smith /* Bottom Level */ 104547c6ae99SBarry Smith for (k=0; k<s_z; k++) { 104647c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 104747c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 1048ce00eea3SSatish Balay x_t = lx[n0 % m]; 104947c6ae99SBarry Smith y_t = ly[(n0 % (m*n))/m]; 105047c6ae99SBarry Smith z_t = lz[n0 / (m*n)]; 105147c6ae99SBarry 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; 10528865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1053ce00eea3SSatish Balay } else if (Xs-xs < 0 && Ys-ys < 0 && Zs-zs < 0) { 10548865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 105547c6ae99SBarry Smith } 105647c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 105747c6ae99SBarry Smith x_t = x; 105847c6ae99SBarry Smith y_t = ly[(n1 % (m*n))/m]; 105947c6ae99SBarry Smith z_t = lz[n1 / (m*n)]; 106047c6ae99SBarry 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; 10618865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1062ce00eea3SSatish Balay } else if (Ys-ys < 0 && Zs-zs < 0) { 10638865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 106447c6ae99SBarry Smith } 106547c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 1066ce00eea3SSatish Balay x_t = lx[n2 % m]; 106747c6ae99SBarry Smith y_t = ly[(n2 % (m*n))/m]; 106847c6ae99SBarry Smith z_t = lz[n2 / (m*n)]; 106947c6ae99SBarry 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; 10708865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1071ce00eea3SSatish Balay } else if (xe-Xe < 0 && Ys-ys < 0 && Zs-zs < 0) { 10728865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 107347c6ae99SBarry Smith } 107447c6ae99SBarry Smith } 107547c6ae99SBarry Smith 107647c6ae99SBarry Smith for (i=0; i<y; i++) { 107747c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 1078ce00eea3SSatish Balay x_t = lx[n3 % m]; 107947c6ae99SBarry Smith y_t = y; 108047c6ae99SBarry Smith z_t = lz[n3 / (m*n)]; 108147c6ae99SBarry Smith s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 10828865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1083ce00eea3SSatish Balay } else if (Xs-xs < 0 && Zs-zs < 0) { 10848865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 108547c6ae99SBarry Smith } 108647c6ae99SBarry Smith 108747c6ae99SBarry Smith if (n4 >= 0) { /* middle */ 108847c6ae99SBarry Smith x_t = x; 108947c6ae99SBarry Smith y_t = y; 109047c6ae99SBarry Smith z_t = lz[n4 / (m*n)]; 109147c6ae99SBarry Smith s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 10928865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1093ce00eea3SSatish Balay } else if (Zs-zs < 0) { 1094bff4a2f0SMatthew G. Knepley if (bz == DM_BOUNDARY_MIRROR) { 1095a6fd6b0aSBarry Smith for (j=0; j<x; j++) idx[nn++] = bases[rank] + j + x*i + (s_z - k - 1)*x*y; 1096c2859e5eSBarry Smith } else { 10978865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 109847c6ae99SBarry Smith } 1099c2859e5eSBarry Smith } 110047c6ae99SBarry Smith 110147c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 1102ce00eea3SSatish Balay x_t = lx[n5 % m]; 110347c6ae99SBarry Smith y_t = y; 110447c6ae99SBarry Smith z_t = lz[n5 / (m*n)]; 110547c6ae99SBarry Smith s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 11068865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1107ce00eea3SSatish Balay } else if (xe-Xe < 0 && Zs-zs < 0) { 11088865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 110947c6ae99SBarry Smith } 111047c6ae99SBarry Smith } 111147c6ae99SBarry Smith 111247c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 111347c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 1114ce00eea3SSatish Balay x_t = lx[n6 % m]; 111547c6ae99SBarry Smith y_t = ly[(n6 % (m*n))/m]; 111647c6ae99SBarry Smith z_t = lz[n6 / (m*n)]; 111747c6ae99SBarry Smith s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 11188865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1119ce00eea3SSatish Balay } else if (Xs-xs < 0 && ye-Ye < 0 && Zs-zs < 0) { 11208865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 112147c6ae99SBarry Smith } 112247c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 112347c6ae99SBarry Smith x_t = x; 112447c6ae99SBarry Smith y_t = ly[(n7 % (m*n))/m]; 112547c6ae99SBarry Smith z_t = lz[n7 / (m*n)]; 112647c6ae99SBarry Smith s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 11278865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1128ce00eea3SSatish Balay } else if (ye-Ye < 0 && Zs-zs < 0) { 11298865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 113047c6ae99SBarry Smith } 113147c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 1132ce00eea3SSatish Balay x_t = lx[n8 % m]; 113347c6ae99SBarry Smith y_t = ly[(n8 % (m*n))/m]; 113447c6ae99SBarry Smith z_t = lz[n8 / (m*n)]; 113547c6ae99SBarry Smith s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 11368865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1137ce00eea3SSatish Balay } else if (xe-Xe < 0 && ye-Ye < 0 && Zs-zs < 0) { 11388865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 113947c6ae99SBarry Smith } 114047c6ae99SBarry Smith } 114147c6ae99SBarry Smith } 114247c6ae99SBarry Smith 114347c6ae99SBarry Smith /* Middle Level */ 114447c6ae99SBarry Smith for (k=0; k<z; k++) { 114547c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 114647c6ae99SBarry Smith if (n9 >= 0) { /* left below */ 1147ce00eea3SSatish Balay x_t = lx[n9 % m]; 114847c6ae99SBarry Smith y_t = ly[(n9 % (m*n))/m]; 114947c6ae99SBarry Smith /* z_t = z; */ 115047c6ae99SBarry Smith s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 11518865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1152ce00eea3SSatish Balay } else if (Xs-xs < 0 && Ys-ys < 0) { 11538865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 115447c6ae99SBarry Smith } 115547c6ae99SBarry Smith if (n10 >= 0) { /* directly below */ 115647c6ae99SBarry Smith x_t = x; 115747c6ae99SBarry Smith y_t = ly[(n10 % (m*n))/m]; 115847c6ae99SBarry Smith /* z_t = z; */ 115947c6ae99SBarry Smith s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 11608865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1161ce00eea3SSatish Balay } else if (Ys-ys < 0) { 1162bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_MIRROR) { 1163feb237baSPierre Jolivet for (j=0; j<x; j++) idx[nn++] = bases[rank] + j + k*x*y + (s_y - i)*x; 1164c2859e5eSBarry Smith } else { 11658865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 1166c2859e5eSBarry Smith } 116747c6ae99SBarry Smith } 116847c6ae99SBarry Smith if (n11 >= 0) { /* right below */ 1169ce00eea3SSatish Balay x_t = lx[n11 % m]; 117047c6ae99SBarry Smith y_t = ly[(n11 % (m*n))/m]; 117147c6ae99SBarry Smith /* z_t = z; */ 117247c6ae99SBarry Smith s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 11738865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1174ce00eea3SSatish Balay } else if (xe-Xe < 0 && Ys-ys < 0) { 11758865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 117647c6ae99SBarry Smith } 117747c6ae99SBarry Smith } 117847c6ae99SBarry Smith 117947c6ae99SBarry Smith for (i=0; i<y; i++) { 118047c6ae99SBarry Smith if (n12 >= 0) { /* directly left */ 1181ce00eea3SSatish Balay x_t = lx[n12 % m]; 118247c6ae99SBarry Smith y_t = y; 118347c6ae99SBarry Smith /* z_t = z; */ 118447c6ae99SBarry Smith s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t; 11858865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1186ce00eea3SSatish Balay } else if (Xs-xs < 0) { 1187bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_MIRROR) { 1188a6fd6b0aSBarry Smith for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + s_x - j - 1 + k*x*y + i*x; 1189c2859e5eSBarry Smith } else { 11908865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 119147c6ae99SBarry Smith } 1192c2859e5eSBarry Smith } 119347c6ae99SBarry Smith 119447c6ae99SBarry Smith /* Interior */ 119547c6ae99SBarry Smith s_t = bases[rank] + i*x + k*x*y; 11968865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = s_t++; 119747c6ae99SBarry Smith 119847c6ae99SBarry Smith if (n14 >= 0) { /* directly right */ 1199ce00eea3SSatish Balay x_t = lx[n14 % m]; 120047c6ae99SBarry Smith y_t = y; 120147c6ae99SBarry Smith /* z_t = z; */ 120247c6ae99SBarry Smith s_t = bases[n14] + i*x_t + k*x_t*y_t; 12038865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1204ce00eea3SSatish Balay } else if (xe-Xe < 0) { 1205bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_MIRROR) { 1206a6fd6b0aSBarry Smith for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x - j - 1 + k*x*y + i*x; 1207c2859e5eSBarry Smith } else { 12088865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 120947c6ae99SBarry Smith } 121047c6ae99SBarry Smith } 1211c2859e5eSBarry Smith } 121247c6ae99SBarry Smith 121347c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 121447c6ae99SBarry Smith if (n15 >= 0) { /* left above */ 1215ce00eea3SSatish Balay x_t = lx[n15 % m]; 121647c6ae99SBarry Smith y_t = ly[(n15 % (m*n))/m]; 121747c6ae99SBarry Smith /* z_t = z; */ 121847c6ae99SBarry Smith s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t; 12198865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1220ce00eea3SSatish Balay } else if (Xs-xs < 0 && ye-Ye < 0) { 12218865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 122247c6ae99SBarry Smith } 122347c6ae99SBarry Smith if (n16 >= 0) { /* directly above */ 122447c6ae99SBarry Smith x_t = x; 122547c6ae99SBarry Smith y_t = ly[(n16 % (m*n))/m]; 122647c6ae99SBarry Smith /* z_t = z; */ 122747c6ae99SBarry Smith s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t; 12288865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1229ce00eea3SSatish Balay } else if (ye-Ye < 0) { 1230bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_MIRROR) { 1231a6fd6b0aSBarry Smith for (j=0; j<x; j++) idx[nn++] = bases[rank] + j + k*x*y + (y-i)*x; 1232c2859e5eSBarry Smith } else { 12338865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 123447c6ae99SBarry Smith } 1235c2859e5eSBarry Smith } 123647c6ae99SBarry Smith if (n17 >= 0) { /* right above */ 1237ce00eea3SSatish Balay x_t = lx[n17 % m]; 123847c6ae99SBarry Smith y_t = ly[(n17 % (m*n))/m]; 123947c6ae99SBarry Smith /* z_t = z; */ 124047c6ae99SBarry Smith s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t; 12418865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1242ce00eea3SSatish Balay } else if (xe-Xe < 0 && ye-Ye < 0) { 12438865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 124447c6ae99SBarry Smith } 124547c6ae99SBarry Smith } 124647c6ae99SBarry Smith } 124747c6ae99SBarry Smith 124847c6ae99SBarry Smith /* Upper Level */ 124947c6ae99SBarry Smith for (k=0; k<s_z; k++) { 125047c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 125147c6ae99SBarry Smith if (n18 >= 0) { /* left below */ 1252ce00eea3SSatish Balay x_t = lx[n18 % m]; 125347c6ae99SBarry Smith y_t = ly[(n18 % (m*n))/m]; 125447c6ae99SBarry Smith /* z_t = lz[n18 / (m*n)]; */ 125547c6ae99SBarry Smith s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 12568865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1257ce00eea3SSatish Balay } else if (Xs-xs < 0 && Ys-ys < 0 && ze-Ze < 0) { 12588865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 125947c6ae99SBarry Smith } 126047c6ae99SBarry Smith if (n19 >= 0) { /* directly below */ 126147c6ae99SBarry Smith x_t = x; 126247c6ae99SBarry Smith y_t = ly[(n19 % (m*n))/m]; 126347c6ae99SBarry Smith /* z_t = lz[n19 / (m*n)]; */ 126447c6ae99SBarry Smith s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 12658865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1266ce00eea3SSatish Balay } else if (Ys-ys < 0 && ze-Ze < 0) { 12678865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 126847c6ae99SBarry Smith } 126947c6ae99SBarry Smith if (n20 >= 0) { /* right below */ 1270ce00eea3SSatish Balay x_t = lx[n20 % m]; 127147c6ae99SBarry Smith y_t = ly[(n20 % (m*n))/m]; 127247c6ae99SBarry Smith /* z_t = lz[n20 / (m*n)]; */ 127347c6ae99SBarry Smith s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 12748865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1275ce00eea3SSatish Balay } else if (xe-Xe < 0 && Ys-ys < 0 && ze-Ze < 0) { 12768865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 127747c6ae99SBarry Smith } 127847c6ae99SBarry Smith } 127947c6ae99SBarry Smith 128047c6ae99SBarry Smith for (i=0; i<y; i++) { 128147c6ae99SBarry Smith if (n21 >= 0) { /* directly left */ 1282ce00eea3SSatish Balay x_t = lx[n21 % m]; 128347c6ae99SBarry Smith y_t = y; 128447c6ae99SBarry Smith /* z_t = lz[n21 / (m*n)]; */ 128547c6ae99SBarry Smith s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t; 12868865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1287ce00eea3SSatish Balay } else if (Xs-xs < 0 && ze-Ze < 0) { 12888865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 128947c6ae99SBarry Smith } 129047c6ae99SBarry Smith 129147c6ae99SBarry Smith if (n22 >= 0) { /* middle */ 129247c6ae99SBarry Smith x_t = x; 129347c6ae99SBarry Smith y_t = y; 129447c6ae99SBarry Smith /* z_t = lz[n22 / (m*n)]; */ 129547c6ae99SBarry Smith s_t = bases[n22] + i*x_t + k*x_t*y_t; 12968865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1297ce00eea3SSatish Balay } else if (ze-Ze < 0) { 1298bff4a2f0SMatthew G. Knepley if (bz == DM_BOUNDARY_MIRROR) { 1299a6fd6b0aSBarry Smith for (j=0; j<x; j++) idx[nn++] = bases[rank] + j + (z-k-1)*x*y + i*x; 1300c2859e5eSBarry Smith } else { 13018865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 130247c6ae99SBarry Smith } 1303c2859e5eSBarry Smith } 130447c6ae99SBarry Smith 130547c6ae99SBarry Smith if (n23 >= 0) { /* directly right */ 1306ce00eea3SSatish Balay x_t = lx[n23 % m]; 130747c6ae99SBarry Smith y_t = y; 130847c6ae99SBarry Smith /* z_t = lz[n23 / (m*n)]; */ 130947c6ae99SBarry Smith s_t = bases[n23] + i*x_t + k*x_t*y_t; 13108865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1311ce00eea3SSatish Balay } else if (xe-Xe < 0 && ze-Ze < 0) { 13128865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 131347c6ae99SBarry Smith } 131447c6ae99SBarry Smith } 131547c6ae99SBarry Smith 131647c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 131747c6ae99SBarry Smith if (n24 >= 0) { /* left above */ 1318ce00eea3SSatish Balay x_t = lx[n24 % m]; 131947c6ae99SBarry Smith y_t = ly[(n24 % (m*n))/m]; 132047c6ae99SBarry Smith /* z_t = lz[n24 / (m*n)]; */ 132147c6ae99SBarry Smith s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t; 13228865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1323ce00eea3SSatish Balay } else if (Xs-xs < 0 && ye-Ye < 0 && ze-Ze < 0) { 13248865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 132547c6ae99SBarry Smith } 132647c6ae99SBarry Smith if (n25 >= 0) { /* directly above */ 132747c6ae99SBarry Smith x_t = x; 132847c6ae99SBarry Smith y_t = ly[(n25 % (m*n))/m]; 132947c6ae99SBarry Smith /* z_t = lz[n25 / (m*n)]; */ 133047c6ae99SBarry Smith s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t; 13318865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1332ce00eea3SSatish Balay } else if (ye-Ye < 0 && ze-Ze < 0) { 13338865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 133447c6ae99SBarry Smith } 133547c6ae99SBarry Smith if (n26 >= 0) { /* right above */ 1336ce00eea3SSatish Balay x_t = lx[n26 % m]; 133747c6ae99SBarry Smith y_t = ly[(n26 % (m*n))/m]; 133847c6ae99SBarry Smith /* z_t = lz[n26 / (m*n)]; */ 133947c6ae99SBarry Smith s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t; 13408865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1341ce00eea3SSatish Balay } else if (xe-Xe < 0 && ye-Ye < 0 && ze-Ze < 0) { 13428865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 134347c6ae99SBarry Smith } 134447c6ae99SBarry Smith } 134547c6ae99SBarry Smith } 134647c6ae99SBarry Smith } 134747c6ae99SBarry Smith /* 134847c6ae99SBarry Smith Set the local to global ordering in the global vector, this allows use 134947c6ae99SBarry Smith of VecSetValuesLocal(). 135047c6ae99SBarry Smith */ 13519566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(comm,dof,nn,idx,PETSC_OWN_POINTER,&da->ltogmap)); 13529566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)da,(PetscObject)da->ltogmap)); 135347c6ae99SBarry Smith 13549566063dSJacob Faibussowitsch PetscCall(PetscFree2(bases,ldims)); 1355ce00eea3SSatish Balay dd->m = m; dd->n = n; dd->p = p; 1356ce00eea3SSatish Balay /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */ 1357ce00eea3SSatish Balay dd->xs = xs*dof; dd->xe = xe*dof; dd->ys = ys; dd->ye = ye; dd->zs = zs; dd->ze = ze; 1358ce00eea3SSatish Balay dd->Xs = Xs*dof; dd->Xe = Xe*dof; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = Zs; dd->Ze = Ze; 1359ce00eea3SSatish Balay 13609566063dSJacob Faibussowitsch PetscCall(VecDestroy(&local)); 13619566063dSJacob Faibussowitsch PetscCall(VecDestroy(&global)); 1362ce00eea3SSatish Balay 1363ce00eea3SSatish Balay dd->gtol = gtol; 1364ce00eea3SSatish Balay dd->base = base; 1365ce00eea3SSatish Balay da->ops->view = DMView_DA_3d; 13660298fd71SBarry Smith dd->ltol = NULL; 13670298fd71SBarry Smith dd->ao = NULL; 136847c6ae99SBarry Smith PetscFunctionReturn(0); 136947c6ae99SBarry Smith } 1370564755cdSBarry Smith 137147c6ae99SBarry Smith /*@C 1372aa219208SBarry Smith DMDACreate3d - Creates an object that will manage the communication of three-dimensional 137347c6ae99SBarry Smith regular array data that is distributed across some processors. 137447c6ae99SBarry Smith 1375d083f849SBarry Smith Collective 137647c6ae99SBarry Smith 137747c6ae99SBarry Smith Input Parameters: 137847c6ae99SBarry Smith + comm - MPI communicator 13791321219cSEthan Coon . bx,by,bz - type of ghost nodes the array have. 1380bff4a2f0SMatthew G. Knepley Use one of DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_PERIODIC. 1381aa219208SBarry Smith . stencil_type - Type of stencil (DMDA_STENCIL_STAR or DMDA_STENCIL_BOX) 1382897f7067SBarry Smith . M,N,P - global dimension in each direction of the array 138347c6ae99SBarry Smith . m,n,p - corresponding number of processors in each dimension 138447c6ae99SBarry Smith (or PETSC_DECIDE to have calculated) 138547c6ae99SBarry Smith . dof - number of degrees of freedom per node 138610d7c030SMatthew G Knepley . s - stencil width 138710d7c030SMatthew G Knepley - lx, ly, lz - arrays containing the number of nodes in each cell along 13880298fd71SBarry Smith the x, y, and z coordinates, or NULL. If non-null, these 138947c6ae99SBarry Smith must be of length as m,n,p and the corresponding 139047c6ae99SBarry Smith m,n, or p cannot be PETSC_DECIDE. Sum of the lx[] entries must be M, sum of 139147c6ae99SBarry Smith the ly[] must N, sum of the lz[] must be P 139247c6ae99SBarry Smith 139347c6ae99SBarry Smith Output Parameter: 139447c6ae99SBarry Smith . da - the resulting distributed array object 139547c6ae99SBarry Smith 139647c6ae99SBarry Smith Options Database Key: 1397706a11cbSBarry Smith + -dm_view - Calls DMView() at the conclusion of DMDACreate3d() 1398e43dc8daSBarry Smith . -da_grid_x <nx> - number of grid points in x direction 1399e43dc8daSBarry Smith . -da_grid_y <ny> - number of grid points in y direction 1400e43dc8daSBarry Smith . -da_grid_z <nz> - number of grid points in z direction 140147c6ae99SBarry Smith . -da_processors_x <MX> - number of processors in x direction 140247c6ae99SBarry Smith . -da_processors_y <MY> - number of processors in y direction 140347c6ae99SBarry Smith . -da_processors_z <MZ> - number of processors in z direction 1404e0f5d30fSBarry Smith . -da_refine_x <rx> - refinement ratio in x direction 1405e0f5d30fSBarry Smith . -da_refine_y <ry> - refinement ratio in y direction 1406e0f5d30fSBarry Smith . -da_refine_z <rz>- refinement ratio in z directio 1407e43dc8daSBarry Smith - -da_refine <n> - refine the DMDA n times before creating it 140847c6ae99SBarry Smith 140947c6ae99SBarry Smith Level: beginner 141047c6ae99SBarry Smith 141147c6ae99SBarry Smith Notes: 1412aa219208SBarry Smith The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the 1413aa219208SBarry Smith standard 7-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes 141447c6ae99SBarry Smith the standard 27-pt stencil. 141547c6ae99SBarry Smith 1416aa219208SBarry Smith The array data itself is NOT stored in the DMDA, it is stored in Vec objects; 1417564755cdSBarry Smith The appropriate vector objects can be obtained with calls to DMCreateGlobalVector() 1418564755cdSBarry Smith and DMCreateLocalVector() and calls to VecDuplicate() if more are needed. 141947c6ae99SBarry Smith 1420897f7067SBarry Smith You must call DMSetUp() after this call before using this DM. 1421897f7067SBarry Smith 1422897f7067SBarry Smith If you wish to use the options database to change values in the DMDA call DMSetFromOptions() after this call 1423897f7067SBarry Smith but before DMSetUp(). 1424897f7067SBarry Smith 1425aa219208SBarry Smith .seealso: DMDestroy(), DMView(), DMDACreate1d(), DMDACreate2d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(), 142699f0b487SRichard Tran Mills DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMLocalToLocalBegin(), DMLocalToLocalEnd(), DMDASetRefinementFactor(), 14273bd220d7SPatrick Sanan DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges(), 14283bd220d7SPatrick Sanan DMStagCreate3d() 142947c6ae99SBarry Smith 143047c6ae99SBarry Smith @*/ 1431bff4a2f0SMatthew G. Knepley PetscErrorCode DMDACreate3d(MPI_Comm comm,DMBoundaryType bx,DMBoundaryType by,DMBoundaryType bz,DMDAStencilType stencil_type,PetscInt M, 14329a42bb27SBarry 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) 143347c6ae99SBarry Smith { 143447c6ae99SBarry Smith PetscFunctionBegin; 14359566063dSJacob Faibussowitsch PetscCall(DMDACreate(comm, da)); 14369566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*da, 3)); 14379566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(*da, M, N, P)); 14389566063dSJacob Faibussowitsch PetscCall(DMDASetNumProcs(*da, m, n, p)); 14399566063dSJacob Faibussowitsch PetscCall(DMDASetBoundaryType(*da, bx, by, bz)); 14409566063dSJacob Faibussowitsch PetscCall(DMDASetDof(*da, dof)); 14419566063dSJacob Faibussowitsch PetscCall(DMDASetStencilType(*da, stencil_type)); 14429566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(*da, s)); 14439566063dSJacob Faibussowitsch PetscCall(DMDASetOwnershipRanges(*da, lx, ly, lz)); 144447c6ae99SBarry Smith PetscFunctionReturn(0); 144547c6ae99SBarry Smith } 1446