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; 15*6858538eSMatthew G. Knepley Vec coordinates; 169a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 179a42bb27SBarry Smith PetscBool ismatlab; 189a42bb27SBarry Smith #endif 1947c6ae99SBarry Smith 2047c6ae99SBarry Smith PetscFunctionBegin; 219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)da),&rank)); 2247c6ae99SBarry Smith 239566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 249566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw)); 259566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERGLVIS,&isglvis)); 269566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary)); 279a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 289566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERMATLAB,&ismatlab)); 299a42bb27SBarry Smith #endif 3047c6ae99SBarry Smith if (iascii) { 3147c6ae99SBarry Smith PetscViewerFormat format; 3247c6ae99SBarry Smith 339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 349566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format)); 3576a8abe0SBarry Smith if (format == PETSC_VIEWER_LOAD_BALANCE) { 3676a8abe0SBarry Smith PetscInt i,nmax = 0,nmin = PETSC_MAX_INT,navg = 0,*nz,nzlocal; 3776a8abe0SBarry Smith DMDALocalInfo info; 3876a8abe0SBarry Smith PetscMPIInt size; 399566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)da),&size)); 409566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da,&info)); 4176a8abe0SBarry Smith nzlocal = info.xm*info.ym*info.zm; 429566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size,&nz)); 439566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&nzlocal,1,MPIU_INT,nz,1,MPIU_INT,PetscObjectComm((PetscObject)da))); 4476a8abe0SBarry Smith for (i=0; i<(PetscInt)size; i++) { 4576a8abe0SBarry Smith nmax = PetscMax(nmax,nz[i]); 4676a8abe0SBarry Smith nmin = PetscMin(nmin,nz[i]); 4776a8abe0SBarry Smith navg += nz[i]; 4876a8abe0SBarry Smith } 499566063dSJacob Faibussowitsch PetscCall(PetscFree(nz)); 5076a8abe0SBarry Smith navg = navg/size; 5163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Load Balance - Grid Points: Min %" PetscInt_FMT " avg %" PetscInt_FMT " max %" PetscInt_FMT "\n",nmin,navg,nmax)); 5276a8abe0SBarry Smith PetscFunctionReturn(0); 5376a8abe0SBarry Smith } 548ec8862eSJed Brown if (format != PETSC_VIEWER_ASCII_VTK_DEPRECATED && format != PETSC_VIEWER_ASCII_VTK_CELL_DEPRECATED && format != PETSC_VIEWER_ASCII_GLVIS) { 55aa219208SBarry Smith DMDALocalInfo info; 569566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da,&info)); 5763a3b9bcSJacob 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)); 5863a3b9bcSJacob 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", 595f80ce2aSJacob Faibussowitsch info.xs,info.xs+info.xm,info.ys,info.ys+info.ym,info.zs,info.zs+info.zm)); 60*6858538eSMatthew G. Knepley PetscCall(DMGetCoordinates(da, &coordinates)); 6147c6ae99SBarry Smith #if !defined(PETSC_USE_COMPLEX) 62*6858538eSMatthew G. Knepley if (coordinates) { 6347c6ae99SBarry Smith PetscInt last; 6447c6ae99SBarry Smith const PetscReal *coors; 65*6858538eSMatthew G. Knepley PetscCall(VecGetArrayRead(coordinates,&coors)); 66*6858538eSMatthew G. Knepley PetscCall(VecGetLocalSize(coordinates,&last)); 6747c6ae99SBarry Smith last = last - 3; 689566063dSJacob 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])); 69*6858538eSMatthew G. Knepley PetscCall(VecRestoreArrayRead(coordinates,&coors)); 7047c6ae99SBarry Smith } 7147c6ae99SBarry Smith #endif 729566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 739566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 741baa6e33SBarry Smith } else if (format == PETSC_VIEWER_ASCII_GLVIS) PetscCall(DMView_DA_GLVis(da,viewer)); 751baa6e33SBarry Smith else PetscCall(DMView_DA_VTK(da,viewer)); 7647c6ae99SBarry Smith } else if (isdraw) { 7747c6ae99SBarry Smith PetscDraw draw; 7847c6ae99SBarry Smith PetscReal ymin = -1.0,ymax = (PetscReal)dd->N; 7947c6ae99SBarry Smith PetscReal xmin = -1.0,xmax = (PetscReal)((dd->M+2)*dd->P),x,y,ycoord,xcoord; 808ea3bf28SBarry Smith PetscInt k,plane,base; 818ea3bf28SBarry Smith const PetscInt *idx; 8247c6ae99SBarry Smith char node[10]; 8347c6ae99SBarry Smith PetscBool isnull; 8447c6ae99SBarry Smith 859566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(viewer,0,&draw)); 869566063dSJacob Faibussowitsch PetscCall(PetscDrawIsNull(draw,&isnull)); 875b399a63SLisandro Dalcin if (isnull) PetscFunctionReturn(0); 8847c6ae99SBarry Smith 899566063dSJacob Faibussowitsch PetscCall(PetscDrawCheckResizedWindow(draw)); 909566063dSJacob Faibussowitsch PetscCall(PetscDrawClear(draw)); 919566063dSJacob Faibussowitsch PetscCall(PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax)); 925b399a63SLisandro Dalcin 93d0609cedSBarry Smith PetscDrawCollectiveBegin(draw); 9447c6ae99SBarry Smith /* first processor draw all node lines */ 95dd400576SPatrick Sanan if (rank == 0) { 9647c6ae99SBarry Smith for (k=0; k<dd->P; k++) { 9747c6ae99SBarry Smith ymin = 0.0; ymax = (PetscReal)(dd->N - 1); 9847c6ae99SBarry Smith for (xmin=(PetscReal)(k*(dd->M+1)); xmin<(PetscReal)(dd->M+(k*(dd->M+1))); xmin++) { 999566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK)); 10047c6ae99SBarry Smith } 10147c6ae99SBarry Smith xmin = (PetscReal)(k*(dd->M+1)); xmax = xmin + (PetscReal)(dd->M - 1); 10247c6ae99SBarry Smith for (ymin=0; ymin<(PetscReal)dd->N; ymin++) { 1039566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK)); 10447c6ae99SBarry Smith } 10547c6ae99SBarry Smith } 10647c6ae99SBarry Smith } 107d0609cedSBarry Smith PetscDrawCollectiveEnd(draw); 1089566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 1099566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw)); 11047c6ae99SBarry Smith 111d0609cedSBarry Smith PetscDrawCollectiveBegin(draw); 1125b399a63SLisandro Dalcin /*Go through and draw for each plane*/ 1135b399a63SLisandro Dalcin for (k=0; k<dd->P; k++) { 11447c6ae99SBarry Smith if ((k >= dd->zs) && (k < dd->ze)) { 11547c6ae99SBarry Smith /* draw my box */ 11647c6ae99SBarry Smith ymin = dd->ys; 11747c6ae99SBarry Smith ymax = dd->ye - 1; 11847c6ae99SBarry Smith xmin = dd->xs/dd->w + (dd->M+1)*k; 11947c6ae99SBarry Smith xmax =(dd->xe-1)/dd->w + (dd->M+1)*k; 12047c6ae99SBarry Smith 1219566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED)); 1229566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED)); 1239566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED)); 1249566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED)); 12547c6ae99SBarry Smith 12647c6ae99SBarry Smith xmin = dd->xs/dd->w; 12747c6ae99SBarry Smith xmax =(dd->xe-1)/dd->w; 12847c6ae99SBarry Smith 129832b7cebSLisandro Dalcin /* identify which processor owns the box */ 1309566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(node,sizeof(node),"%d",(int)rank)); 1319566063dSJacob Faibussowitsch PetscCall(PetscDrawString(draw,xmin+(dd->M+1)*k+.2,ymin+.3,PETSC_DRAW_RED,node)); 13247c6ae99SBarry Smith /* put in numbers*/ 13347c6ae99SBarry Smith base = (dd->base+(dd->xe-dd->xs)*(dd->ye-dd->ys)*(k-dd->zs))/dd->w; 13447c6ae99SBarry Smith for (y=ymin; y<=ymax; y++) { 13547c6ae99SBarry Smith for (x=xmin+(dd->M+1)*k; x<=xmax+(dd->M+1)*k; x++) { 1369566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(node,sizeof(node),"%d",(int)base++)); 1379566063dSJacob Faibussowitsch PetscCall(PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node)); 13847c6ae99SBarry Smith } 13947c6ae99SBarry Smith } 14047c6ae99SBarry Smith 14147c6ae99SBarry Smith } 14247c6ae99SBarry Smith } 143d0609cedSBarry Smith PetscDrawCollectiveEnd(draw); 1449566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 1459566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw)); 14647c6ae99SBarry Smith 147d0609cedSBarry Smith PetscDrawCollectiveBegin(draw); 14847c6ae99SBarry Smith for (k=0-dd->s; k<dd->P+dd->s; k++) { 14947c6ae99SBarry Smith /* Go through and draw for each plane */ 15047c6ae99SBarry Smith if ((k >= dd->Zs) && (k < dd->Ze)) { 15147c6ae99SBarry Smith /* overlay ghost numbers, useful for error checking */ 152565245c5SBarry Smith base = (dd->Xe-dd->Xs)*(dd->Ye-dd->Ys)*(k-dd->Zs)/dd->w; 1539566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(da->ltogmap,&idx)); 15447c6ae99SBarry Smith plane=k; 155565245c5SBarry Smith /* Keep z wrap around points on the drawing */ 1568865f1eaSKarl Rupp if (k<0) plane=dd->P+k; 1578865f1eaSKarl Rupp if (k>=dd->P) plane=k-dd->P; 15847c6ae99SBarry Smith ymin = dd->Ys; ymax = dd->Ye; 15947c6ae99SBarry Smith xmin = (dd->M+1)*plane*dd->w; 16047c6ae99SBarry Smith xmax = (dd->M+1)*plane*dd->w+dd->M*dd->w; 16147c6ae99SBarry Smith for (y=ymin; y<ymax; y++) { 16247c6ae99SBarry Smith for (x=xmin+dd->Xs; x<xmin+dd->Xe; x+=dd->w) { 1638ea3bf28SBarry Smith sprintf(node,"%d",(int)(idx[base])); 16447c6ae99SBarry Smith ycoord = y; 16547c6ae99SBarry Smith /*Keep y wrap around points on drawing */ 1668865f1eaSKarl Rupp if (y<0) ycoord = dd->N+y; 1678865f1eaSKarl Rupp if (y>=dd->N) ycoord = y-dd->N; 16847c6ae99SBarry Smith xcoord = x; /* Keep x wrap points on drawing */ 1698865f1eaSKarl Rupp if (x<xmin) xcoord = xmax - (xmin-x); 1708865f1eaSKarl Rupp if (x>=xmax) xcoord = xmin + (x-xmax); 1719566063dSJacob Faibussowitsch PetscCall(PetscDrawString(draw,xcoord/dd->w,ycoord,PETSC_DRAW_BLUE,node)); 172565245c5SBarry Smith base++; 17347c6ae99SBarry Smith } 17447c6ae99SBarry Smith } 1759566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(da->ltogmap,&idx)); 17647c6ae99SBarry Smith } 17747c6ae99SBarry Smith } 178d0609cedSBarry Smith PetscDrawCollectiveEnd(draw); 1799566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 1809566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw)); 1819566063dSJacob Faibussowitsch PetscCall(PetscDrawSave(draw)); 182c9493c35SStefano Zampini } else if (isglvis) { 1839566063dSJacob Faibussowitsch PetscCall(DMView_DA_GLVis(da,viewer)); 1849a42bb27SBarry Smith } else if (isbinary) { 1859566063dSJacob Faibussowitsch PetscCall(DMView_DA_Binary(da,viewer)); 1869a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 1879a42bb27SBarry Smith } else if (ismatlab) { 1889566063dSJacob Faibussowitsch PetscCall(DMView_DA_Matlab(da,viewer)); 1899a42bb27SBarry Smith #endif 19011aeaf0aSBarry Smith } 19147c6ae99SBarry Smith PetscFunctionReturn(0); 19247c6ae99SBarry Smith } 19347c6ae99SBarry Smith 1947087cfbeSBarry Smith PetscErrorCode DMSetUp_DA_3D(DM da) 19547c6ae99SBarry Smith { 19647c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 19747c6ae99SBarry Smith const PetscInt M = dd->M; 19847c6ae99SBarry Smith const PetscInt N = dd->N; 19947c6ae99SBarry Smith const PetscInt P = dd->P; 20047c6ae99SBarry Smith PetscInt m = dd->m; 20147c6ae99SBarry Smith PetscInt n = dd->n; 20247c6ae99SBarry Smith PetscInt p = dd->p; 20347c6ae99SBarry Smith const PetscInt dof = dd->w; 20447c6ae99SBarry Smith const PetscInt s = dd->s; 205bff4a2f0SMatthew G. Knepley DMBoundaryType bx = dd->bx; 206bff4a2f0SMatthew G. Knepley DMBoundaryType by = dd->by; 207bff4a2f0SMatthew G. Knepley DMBoundaryType bz = dd->bz; 20819fd82e9SBarry Smith DMDAStencilType stencil_type = dd->stencil_type; 20947c6ae99SBarry Smith PetscInt *lx = dd->lx; 21047c6ae99SBarry Smith PetscInt *ly = dd->ly; 21147c6ae99SBarry Smith PetscInt *lz = dd->lz; 21247c6ae99SBarry Smith MPI_Comm comm; 21347c6ae99SBarry Smith PetscMPIInt rank,size; 214ce00eea3SSatish Balay PetscInt xs = 0,xe,ys = 0,ye,zs = 0,ze,x = 0,y = 0,z = 0; 215bd1fc5aeSBarry Smith PetscInt Xs,Xe,Ys,Ye,Zs,Ze,IXs,IXe,IYs,IYe,IZs,IZe,pm; 2168ea3bf28SBarry Smith PetscInt left,right,up,down,bottom,top,i,j,k,*idx,nn; 21747c6ae99SBarry Smith PetscInt n0,n1,n2,n3,n4,n5,n6,n7,n8,n9,n10,n11,n12,n14; 21847c6ae99SBarry Smith PetscInt n15,n16,n17,n18,n19,n20,n21,n22,n23,n24,n25,n26; 219db87c5ecSEthan Coon PetscInt *bases,*ldims,base,x_t,y_t,z_t,s_t,count,s_x,s_y,s_z; 22047c6ae99SBarry Smith PetscInt sn0 = 0,sn1 = 0,sn2 = 0,sn3 = 0,sn5 = 0,sn6 = 0,sn7 = 0; 22147c6ae99SBarry Smith PetscInt sn8 = 0,sn9 = 0,sn11 = 0,sn15 = 0,sn24 = 0,sn25 = 0,sn26 = 0; 22247c6ae99SBarry Smith PetscInt sn17 = 0,sn18 = 0,sn19 = 0,sn20 = 0,sn21 = 0,sn23 = 0; 22347c6ae99SBarry Smith Vec local,global; 224bd1fc5aeSBarry Smith VecScatter gtol; 22545b6f7e9SBarry Smith IS to,from; 2266f951b95Secoon PetscBool twod; 22747c6ae99SBarry Smith 22847c6ae99SBarry Smith PetscFunctionBegin; 2297a8be351SBarry 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"); 2309566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject) da, &comm)); 2313855c12bSBarry Smith #if !defined(PETSC_USE_64BIT_INDICES) 23263a3b9bcSJacob 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); 2333855c12bSBarry Smith #endif 2343855c12bSBarry Smith 2359566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 2369566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm,&rank)); 23747c6ae99SBarry Smith 23847c6ae99SBarry Smith if (m != PETSC_DECIDE) { 23963a3b9bcSJacob Faibussowitsch PetscCheck(m >= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %" PetscInt_FMT,m); 24063a3b9bcSJacob Faibussowitsch else PetscCheck(m <= size,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %" PetscInt_FMT " %d",m,size); 24147c6ae99SBarry Smith } 24247c6ae99SBarry Smith if (n != PETSC_DECIDE) { 24363a3b9bcSJacob Faibussowitsch PetscCheck(n >= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %" PetscInt_FMT,n); 24463a3b9bcSJacob Faibussowitsch else PetscCheck(n <= size,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %" PetscInt_FMT " %d",n,size); 24547c6ae99SBarry Smith } 24647c6ae99SBarry Smith if (p != PETSC_DECIDE) { 24763a3b9bcSJacob Faibussowitsch PetscCheck(p >= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Z direction: %" PetscInt_FMT,p); 24863a3b9bcSJacob Faibussowitsch else PetscCheck(p <= size,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Z direction: %" PetscInt_FMT " %d",p,size); 24947c6ae99SBarry Smith } 2501dca8a05SBarry 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); 25147c6ae99SBarry Smith 25247c6ae99SBarry Smith /* Partition the array among the processors */ 25347c6ae99SBarry Smith if (m == PETSC_DECIDE && n != PETSC_DECIDE && p != PETSC_DECIDE) { 25447c6ae99SBarry Smith m = size/(n*p); 25547c6ae99SBarry Smith } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) { 25647c6ae99SBarry Smith n = size/(m*p); 25747c6ae99SBarry Smith } else if (m != PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) { 25847c6ae99SBarry Smith p = size/(m*n); 25947c6ae99SBarry Smith } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) { 26047c6ae99SBarry Smith /* try for squarish distribution */ 261369cc0aeSBarry Smith m = (int)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)N*p))); 26247c6ae99SBarry Smith if (!m) m = 1; 26347c6ae99SBarry Smith while (m > 0) { 26447c6ae99SBarry Smith n = size/(m*p); 26547c6ae99SBarry Smith if (m*n*p == size) break; 26647c6ae99SBarry Smith m--; 26747c6ae99SBarry Smith } 26863a3b9bcSJacob Faibussowitsch PetscCheck(m,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad p value: p = %" PetscInt_FMT,p); 26947c6ae99SBarry Smith if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;} 27047c6ae99SBarry Smith } else if (m == PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) { 27147c6ae99SBarry Smith /* try for squarish distribution */ 272369cc0aeSBarry Smith m = (int)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)P*n))); 27347c6ae99SBarry Smith if (!m) m = 1; 27447c6ae99SBarry Smith while (m > 0) { 27547c6ae99SBarry Smith p = size/(m*n); 27647c6ae99SBarry Smith if (m*n*p == size) break; 27747c6ae99SBarry Smith m--; 27847c6ae99SBarry Smith } 27963a3b9bcSJacob Faibussowitsch PetscCheck(m,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad n value: n = %" PetscInt_FMT,n); 28047c6ae99SBarry Smith if (M > P && m < p) {PetscInt _m = m; m = p; p = _m;} 28147c6ae99SBarry Smith } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) { 28247c6ae99SBarry Smith /* try for squarish distribution */ 283369cc0aeSBarry Smith n = (int)(0.5 + PetscSqrtReal(((PetscReal)N)*((PetscReal)size)/((PetscReal)P*m))); 28447c6ae99SBarry Smith if (!n) n = 1; 28547c6ae99SBarry Smith while (n > 0) { 28647c6ae99SBarry Smith p = size/(m*n); 28747c6ae99SBarry Smith if (m*n*p == size) break; 28847c6ae99SBarry Smith n--; 28947c6ae99SBarry Smith } 29063a3b9bcSJacob Faibussowitsch PetscCheck(n,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad m value: m = %" PetscInt_FMT,n); 29147c6ae99SBarry Smith if (N > P && n < p) {PetscInt _n = n; n = p; p = _n;} 29247c6ae99SBarry Smith } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) { 29347c6ae99SBarry Smith /* try for squarish distribution */ 294369cc0aeSBarry Smith n = (PetscInt)(0.5 + PetscPowReal(((PetscReal)N*N)*((PetscReal)size)/((PetscReal)P*M),(PetscReal)(1./3.))); 29547c6ae99SBarry Smith if (!n) n = 1; 29647c6ae99SBarry Smith while (n > 0) { 29747c6ae99SBarry Smith pm = size/n; 29847c6ae99SBarry Smith if (n*pm == size) break; 29947c6ae99SBarry Smith n--; 30047c6ae99SBarry Smith } 30147c6ae99SBarry Smith if (!n) n = 1; 302369cc0aeSBarry Smith m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)P*n))); 30347c6ae99SBarry Smith if (!m) m = 1; 30447c6ae99SBarry Smith while (m > 0) { 30547c6ae99SBarry Smith p = size/(m*n); 30647c6ae99SBarry Smith if (m*n*p == size) break; 30747c6ae99SBarry Smith m--; 30847c6ae99SBarry Smith } 30947c6ae99SBarry Smith if (M > P && m < p) {PetscInt _m = m; m = p; p = _m;} 31008401ef6SPierre Jolivet } else PetscCheck(m*n*p == size,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition"); 31147c6ae99SBarry Smith 31208401ef6SPierre Jolivet PetscCheck(m*n*p == size,PetscObjectComm((PetscObject)da),PETSC_ERR_PLIB,"Could not find good partition"); 31363a3b9bcSJacob Faibussowitsch PetscCheck(M >= m,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %" PetscInt_FMT " %" PetscInt_FMT,M,m); 31463a3b9bcSJacob Faibussowitsch PetscCheck(N >= n,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %" PetscInt_FMT " %" PetscInt_FMT,N,n); 31563a3b9bcSJacob Faibussowitsch PetscCheck(P >= p,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Partition in z direction is too fine! %" PetscInt_FMT " %" PetscInt_FMT,P,p); 31647c6ae99SBarry Smith 31747c6ae99SBarry Smith /* 31847c6ae99SBarry Smith Determine locally owned region 31947c6ae99SBarry Smith [x, y, or z]s is the first local node number, [x, y, z] is the number of local nodes 32047c6ae99SBarry Smith */ 32147c6ae99SBarry Smith 32247c6ae99SBarry Smith if (!lx) { 3239566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m, &dd->lx)); 32447c6ae99SBarry Smith lx = dd->lx; 3258865f1eaSKarl Rupp for (i=0; i<m; i++) lx[i] = M/m + ((M % m) > (i % m)); 32647c6ae99SBarry Smith } 32747c6ae99SBarry Smith x = lx[rank % m]; 32847c6ae99SBarry Smith xs = 0; 3298865f1eaSKarl Rupp for (i=0; i<(rank%m); i++) xs += lx[i]; 3301dca8a05SBarry 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); 33147c6ae99SBarry Smith 33247c6ae99SBarry Smith if (!ly) { 3339566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &dd->ly)); 33447c6ae99SBarry Smith ly = dd->ly; 3358865f1eaSKarl Rupp for (i=0; i<n; i++) ly[i] = N/n + ((N % n) > (i % n)); 33647c6ae99SBarry Smith } 33747c6ae99SBarry Smith y = ly[(rank % (m*n))/m]; 3381dca8a05SBarry 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); 33930729d88SBarry Smith 34047c6ae99SBarry Smith ys = 0; 3418865f1eaSKarl Rupp for (i=0; i<(rank % (m*n))/m; i++) ys += ly[i]; 34247c6ae99SBarry Smith 34347c6ae99SBarry Smith if (!lz) { 3449566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(p, &dd->lz)); 34547c6ae99SBarry Smith lz = dd->lz; 3468865f1eaSKarl Rupp for (i=0; i<p; i++) lz[i] = P/p + ((P % p) > (i % p)); 34747c6ae99SBarry Smith } 34847c6ae99SBarry Smith z = lz[rank/(m*n)]; 349bcea557cSEthan Coon 350fdc81ce8SEthan Coon /* note this is different than x- and y-, as we will handle as an important special 351bff4a2f0SMatthew G. Knepley case when p=P=1 and DM_BOUNDARY_PERIODIC and s > z. This is to deal with 2D problems 352fdc81ce8SEthan Coon in a 3D code. Additional code for this case is noted with "2d case" comments */ 3536f951b95Secoon twod = PETSC_FALSE; 3548865f1eaSKarl Rupp if (P == 1) twod = PETSC_TRUE; 3551dca8a05SBarry 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); 35647c6ae99SBarry Smith zs = 0; 3578865f1eaSKarl Rupp for (i=0; i<(rank/(m*n)); i++) zs += lz[i]; 35847c6ae99SBarry Smith ye = ys + y; 35947c6ae99SBarry Smith xe = xs + x; 36047c6ae99SBarry Smith ze = zs + z; 36147c6ae99SBarry Smith 36288661749SPeter Brune /* determine ghost region (Xs) and region scattered into (IXs) */ 363d9c9ebe5SPeter Brune if (xs-s > 0) { 364d9c9ebe5SPeter Brune Xs = xs - s; IXs = xs - s; 36588661749SPeter Brune } else { 3668865f1eaSKarl Rupp if (bx) Xs = xs - s; 3678865f1eaSKarl Rupp else Xs = 0; 36888661749SPeter Brune IXs = 0; 36988661749SPeter Brune } 370d9c9ebe5SPeter Brune if (xe+s <= M) { 371d9c9ebe5SPeter Brune Xe = xe + s; IXe = xe + s; 37288661749SPeter Brune } else { 37388661749SPeter Brune if (bx) { 374d9c9ebe5SPeter Brune Xs = xs - s; Xe = xe + s; 3758865f1eaSKarl Rupp } else Xe = M; 37688661749SPeter Brune IXe = M; 37788661749SPeter Brune } 37847c6ae99SBarry Smith 379bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC || bx == DM_BOUNDARY_MIRROR) { 380d9c9ebe5SPeter Brune IXs = xs - s; 381d9c9ebe5SPeter Brune IXe = xe + s; 382d9c9ebe5SPeter Brune Xs = xs - s; 383d9c9ebe5SPeter Brune Xe = xe + s; 38488661749SPeter Brune } 38588661749SPeter Brune 386d9c9ebe5SPeter Brune if (ys-s > 0) { 387d9c9ebe5SPeter Brune Ys = ys - s; IYs = ys - s; 38888661749SPeter Brune } else { 3898865f1eaSKarl Rupp if (by) Ys = ys - s; 3908865f1eaSKarl Rupp else Ys = 0; 39188661749SPeter Brune IYs = 0; 39288661749SPeter Brune } 393d9c9ebe5SPeter Brune if (ye+s <= N) { 394d9c9ebe5SPeter Brune Ye = ye + s; IYe = ye + s; 39588661749SPeter Brune } else { 3968865f1eaSKarl Rupp if (by) Ye = ye + s; 3978865f1eaSKarl Rupp else Ye = N; 39888661749SPeter Brune IYe = N; 39988661749SPeter Brune } 40088661749SPeter Brune 401bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_PERIODIC || by == DM_BOUNDARY_MIRROR) { 402d9c9ebe5SPeter Brune IYs = ys - s; 403d9c9ebe5SPeter Brune IYe = ye + s; 404d9c9ebe5SPeter Brune Ys = ys - s; 405d9c9ebe5SPeter Brune Ye = ye + s; 40688661749SPeter Brune } 40788661749SPeter Brune 408d9c9ebe5SPeter Brune if (zs-s > 0) { 409d9c9ebe5SPeter Brune Zs = zs - s; IZs = zs - s; 41088661749SPeter Brune } else { 4118865f1eaSKarl Rupp if (bz) Zs = zs - s; 4128865f1eaSKarl Rupp else Zs = 0; 41388661749SPeter Brune IZs = 0; 41488661749SPeter Brune } 415d9c9ebe5SPeter Brune if (ze+s <= P) { 416d9c9ebe5SPeter Brune Ze = ze + s; IZe = ze + s; 41788661749SPeter Brune } else { 4188865f1eaSKarl Rupp if (bz) Ze = ze + s; 4198865f1eaSKarl Rupp else Ze = P; 42088661749SPeter Brune IZe = P; 42188661749SPeter Brune } 42288661749SPeter Brune 423bff4a2f0SMatthew G. Knepley if (bz == DM_BOUNDARY_PERIODIC || bz == DM_BOUNDARY_MIRROR) { 424d9c9ebe5SPeter Brune IZs = zs - s; 425d9c9ebe5SPeter Brune IZe = ze + s; 426d9c9ebe5SPeter Brune Zs = zs - s; 427d9c9ebe5SPeter Brune Ze = ze + s; 42888661749SPeter Brune } 42947c6ae99SBarry Smith 43047c6ae99SBarry Smith /* Resize all X parameters to reflect w */ 431d9c9ebe5SPeter Brune s_x = s; 432d9c9ebe5SPeter Brune s_y = s; 433d9c9ebe5SPeter Brune s_z = s; 43447c6ae99SBarry Smith 43547c6ae99SBarry Smith /* determine starting point of each processor */ 43647c6ae99SBarry Smith nn = x*y*z; 4379566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(size+1,&bases,size,&ldims)); 4389566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm)); 43947c6ae99SBarry Smith bases[0] = 0; 4408865f1eaSKarl Rupp for (i=1; i<=size; i++) bases[i] = ldims[i-1]; 4418865f1eaSKarl Rupp for (i=1; i<=size; i++) bases[i] += bases[i-1]; 442ce00eea3SSatish Balay base = bases[rank]*dof; 44347c6ae99SBarry Smith 44447c6ae99SBarry Smith /* allocate the base parallel and sequential vectors */ 445ce00eea3SSatish Balay dd->Nlocal = x*y*z*dof; 4469566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,NULL,&global)); 447ce00eea3SSatish Balay dd->nlocal = (Xe-Xs)*(Ye-Ys)*(Ze-Zs)*dof; 4489566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,NULL,&local)); 44947c6ae99SBarry Smith 4504104a7a0SPatrick Sanan /* generate global to local vector scatter and local to global mapping*/ 45147c6ae99SBarry Smith 452ce00eea3SSatish Balay /* global to local must include ghost points within the domain, 453ce00eea3SSatish Balay but not ghost points outside the domain that aren't periodic */ 4549566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((IXe-IXs)*(IYe-IYs)*(IZe-IZs),&idx)); 455d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_BOX) { 456ce00eea3SSatish Balay left = IXs - Xs; right = left + (IXe-IXs); 457ce00eea3SSatish Balay bottom = IYs - Ys; top = bottom + (IYe-IYs); 458ce00eea3SSatish Balay down = IZs - Zs; up = down + (IZe-IZs); 459ce00eea3SSatish Balay count = 0; 460ce00eea3SSatish Balay for (i=down; i<up; i++) { 461ce00eea3SSatish Balay for (j=bottom; j<top; j++) { 462ce00eea3SSatish Balay for (k=left; k<right; k++) { 463ce00eea3SSatish Balay idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 464ce00eea3SSatish Balay } 465ce00eea3SSatish Balay } 466ce00eea3SSatish Balay } 4679566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to)); 46847c6ae99SBarry Smith } else { 46947c6ae99SBarry Smith /* This is way ugly! We need to list the funny cross type region */ 470ce00eea3SSatish Balay left = xs - Xs; right = left + x; 47147c6ae99SBarry Smith bottom = ys - Ys; top = bottom + y; 47247c6ae99SBarry Smith down = zs - Zs; up = down + z; 47347c6ae99SBarry Smith count = 0; 474a5b23f4aSJose E. Roman /* the bottom chunk */ 475ce00eea3SSatish Balay for (i=(IZs-Zs); i<down; i++) { 47647c6ae99SBarry Smith for (j=bottom; j<top; j++) { 477ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 47847c6ae99SBarry Smith } 47947c6ae99SBarry Smith } 48047c6ae99SBarry Smith /* the middle piece */ 48147c6ae99SBarry Smith for (i=down; i<up; i++) { 48247c6ae99SBarry Smith /* front */ 483ce00eea3SSatish Balay for (j=(IYs-Ys); j<bottom; j++) { 484ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 48547c6ae99SBarry Smith } 48647c6ae99SBarry Smith /* middle */ 48747c6ae99SBarry Smith for (j=bottom; j<top; j++) { 488ce00eea3SSatish Balay for (k=IXs-Xs; k<IXe-Xs; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 48947c6ae99SBarry Smith } 49047c6ae99SBarry Smith /* back */ 491ce00eea3SSatish Balay for (j=top; j<top+IYe-ye; j++) { 492ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 49347c6ae99SBarry Smith } 49447c6ae99SBarry Smith } 49547c6ae99SBarry Smith /* the top piece */ 496ce00eea3SSatish Balay for (i=up; i<up+IZe-ze; i++) { 49747c6ae99SBarry Smith for (j=bottom; j<top; j++) { 498ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 49947c6ae99SBarry Smith } 50047c6ae99SBarry Smith } 5019566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to)); 50247c6ae99SBarry Smith } 50347c6ae99SBarry Smith 50447c6ae99SBarry Smith /* determine who lies on each side of use stored in n24 n25 n26 50547c6ae99SBarry Smith n21 n22 n23 50647c6ae99SBarry Smith n18 n19 n20 50747c6ae99SBarry Smith 50847c6ae99SBarry Smith n15 n16 n17 50947c6ae99SBarry Smith n12 n14 51047c6ae99SBarry Smith n9 n10 n11 51147c6ae99SBarry Smith 51247c6ae99SBarry Smith n6 n7 n8 51347c6ae99SBarry Smith n3 n4 n5 51447c6ae99SBarry Smith n0 n1 n2 51547c6ae99SBarry Smith */ 51647c6ae99SBarry Smith 51747c6ae99SBarry Smith /* Solve for X,Y, and Z Periodic Case First, Then Modify Solution */ 51847c6ae99SBarry Smith /* Assume Nodes are Internal to the Cube */ 51947c6ae99SBarry Smith n0 = rank - m*n - m - 1; 52047c6ae99SBarry Smith n1 = rank - m*n - m; 52147c6ae99SBarry Smith n2 = rank - m*n - m + 1; 52247c6ae99SBarry Smith n3 = rank - m*n -1; 52347c6ae99SBarry Smith n4 = rank - m*n; 52447c6ae99SBarry Smith n5 = rank - m*n + 1; 52547c6ae99SBarry Smith n6 = rank - m*n + m - 1; 52647c6ae99SBarry Smith n7 = rank - m*n + m; 52747c6ae99SBarry Smith n8 = rank - m*n + m + 1; 52847c6ae99SBarry Smith 52947c6ae99SBarry Smith n9 = rank - m - 1; 53047c6ae99SBarry Smith n10 = rank - m; 53147c6ae99SBarry Smith n11 = rank - m + 1; 53247c6ae99SBarry Smith n12 = rank - 1; 53347c6ae99SBarry Smith n14 = rank + 1; 53447c6ae99SBarry Smith n15 = rank + m - 1; 53547c6ae99SBarry Smith n16 = rank + m; 53647c6ae99SBarry Smith n17 = rank + m + 1; 53747c6ae99SBarry Smith 53847c6ae99SBarry Smith n18 = rank + m*n - m - 1; 53947c6ae99SBarry Smith n19 = rank + m*n - m; 54047c6ae99SBarry Smith n20 = rank + m*n - m + 1; 54147c6ae99SBarry Smith n21 = rank + m*n - 1; 54247c6ae99SBarry Smith n22 = rank + m*n; 54347c6ae99SBarry Smith n23 = rank + m*n + 1; 54447c6ae99SBarry Smith n24 = rank + m*n + m - 1; 54547c6ae99SBarry Smith n25 = rank + m*n + m; 54647c6ae99SBarry Smith n26 = rank + m*n + m + 1; 54747c6ae99SBarry Smith 54847c6ae99SBarry Smith /* Assume Pieces are on Faces of Cube */ 54947c6ae99SBarry Smith 55047c6ae99SBarry Smith if (xs == 0) { /* First assume not corner or edge */ 55147c6ae99SBarry Smith n0 = rank -1 - (m*n); 55247c6ae99SBarry Smith n3 = rank + m -1 - (m*n); 55347c6ae99SBarry Smith n6 = rank + 2*m -1 - (m*n); 55447c6ae99SBarry Smith n9 = rank -1; 55547c6ae99SBarry Smith n12 = rank + m -1; 55647c6ae99SBarry Smith n15 = rank + 2*m -1; 55747c6ae99SBarry Smith n18 = rank -1 + (m*n); 55847c6ae99SBarry Smith n21 = rank + m -1 + (m*n); 55947c6ae99SBarry Smith n24 = rank + 2*m -1 + (m*n); 56047c6ae99SBarry Smith } 56147c6ae99SBarry Smith 562ce00eea3SSatish Balay if (xe == M) { /* First assume not corner or edge */ 56347c6ae99SBarry Smith n2 = rank -2*m +1 - (m*n); 56447c6ae99SBarry Smith n5 = rank - m +1 - (m*n); 56547c6ae99SBarry Smith n8 = rank +1 - (m*n); 56647c6ae99SBarry Smith n11 = rank -2*m +1; 56747c6ae99SBarry Smith n14 = rank - m +1; 56847c6ae99SBarry Smith n17 = rank +1; 56947c6ae99SBarry Smith n20 = rank -2*m +1 + (m*n); 57047c6ae99SBarry Smith n23 = rank - m +1 + (m*n); 57147c6ae99SBarry Smith n26 = rank +1 + (m*n); 57247c6ae99SBarry Smith } 57347c6ae99SBarry Smith 57447c6ae99SBarry Smith if (ys==0) { /* First assume not corner or edge */ 57547c6ae99SBarry Smith n0 = rank + m * (n-1) -1 - (m*n); 57647c6ae99SBarry Smith n1 = rank + m * (n-1) - (m*n); 57747c6ae99SBarry Smith n2 = rank + m * (n-1) +1 - (m*n); 57847c6ae99SBarry Smith n9 = rank + m * (n-1) -1; 57947c6ae99SBarry Smith n10 = rank + m * (n-1); 58047c6ae99SBarry Smith n11 = rank + m * (n-1) +1; 58147c6ae99SBarry Smith n18 = rank + m * (n-1) -1 + (m*n); 58247c6ae99SBarry Smith n19 = rank + m * (n-1) + (m*n); 58347c6ae99SBarry Smith n20 = rank + m * (n-1) +1 + (m*n); 58447c6ae99SBarry Smith } 58547c6ae99SBarry Smith 58647c6ae99SBarry Smith if (ye == N) { /* First assume not corner or edge */ 58747c6ae99SBarry Smith n6 = rank - m * (n-1) -1 - (m*n); 58847c6ae99SBarry Smith n7 = rank - m * (n-1) - (m*n); 58947c6ae99SBarry Smith n8 = rank - m * (n-1) +1 - (m*n); 59047c6ae99SBarry Smith n15 = rank - m * (n-1) -1; 59147c6ae99SBarry Smith n16 = rank - m * (n-1); 59247c6ae99SBarry Smith n17 = rank - m * (n-1) +1; 59347c6ae99SBarry Smith n24 = rank - m * (n-1) -1 + (m*n); 59447c6ae99SBarry Smith n25 = rank - m * (n-1) + (m*n); 59547c6ae99SBarry Smith n26 = rank - m * (n-1) +1 + (m*n); 59647c6ae99SBarry Smith } 59747c6ae99SBarry Smith 59847c6ae99SBarry Smith if (zs == 0) { /* First assume not corner or edge */ 59947c6ae99SBarry Smith n0 = size - (m*n) + rank - m - 1; 60047c6ae99SBarry Smith n1 = size - (m*n) + rank - m; 60147c6ae99SBarry Smith n2 = size - (m*n) + rank - m + 1; 60247c6ae99SBarry Smith n3 = size - (m*n) + rank - 1; 60347c6ae99SBarry Smith n4 = size - (m*n) + rank; 60447c6ae99SBarry Smith n5 = size - (m*n) + rank + 1; 60547c6ae99SBarry Smith n6 = size - (m*n) + rank + m - 1; 60647c6ae99SBarry Smith n7 = size - (m*n) + rank + m; 60747c6ae99SBarry Smith n8 = size - (m*n) + rank + m + 1; 60847c6ae99SBarry Smith } 60947c6ae99SBarry Smith 61047c6ae99SBarry Smith if (ze == P) { /* First assume not corner or edge */ 61147c6ae99SBarry Smith n18 = (m*n) - (size-rank) - m - 1; 61247c6ae99SBarry Smith n19 = (m*n) - (size-rank) - m; 61347c6ae99SBarry Smith n20 = (m*n) - (size-rank) - m + 1; 61447c6ae99SBarry Smith n21 = (m*n) - (size-rank) - 1; 61547c6ae99SBarry Smith n22 = (m*n) - (size-rank); 61647c6ae99SBarry Smith n23 = (m*n) - (size-rank) + 1; 61747c6ae99SBarry Smith n24 = (m*n) - (size-rank) + m - 1; 61847c6ae99SBarry Smith n25 = (m*n) - (size-rank) + m; 61947c6ae99SBarry Smith n26 = (m*n) - (size-rank) + m + 1; 62047c6ae99SBarry Smith } 62147c6ae99SBarry Smith 62247c6ae99SBarry Smith if ((xs==0) && (zs==0)) { /* Assume an edge, not corner */ 62347c6ae99SBarry Smith n0 = size - m*n + rank + m-1 - m; 62447c6ae99SBarry Smith n3 = size - m*n + rank + m-1; 62547c6ae99SBarry Smith n6 = size - m*n + rank + m-1 + m; 62647c6ae99SBarry Smith } 62747c6ae99SBarry Smith 62847c6ae99SBarry Smith if ((xs==0) && (ze==P)) { /* Assume an edge, not corner */ 62947c6ae99SBarry Smith n18 = m*n - (size - rank) + m-1 - m; 63047c6ae99SBarry Smith n21 = m*n - (size - rank) + m-1; 63147c6ae99SBarry Smith n24 = m*n - (size - rank) + m-1 + m; 63247c6ae99SBarry Smith } 63347c6ae99SBarry Smith 63447c6ae99SBarry Smith if ((xs==0) && (ys==0)) { /* Assume an edge, not corner */ 63547c6ae99SBarry Smith n0 = rank + m*n -1 - m*n; 63647c6ae99SBarry Smith n9 = rank + m*n -1; 63747c6ae99SBarry Smith n18 = rank + m*n -1 + m*n; 63847c6ae99SBarry Smith } 63947c6ae99SBarry Smith 64047c6ae99SBarry Smith if ((xs==0) && (ye==N)) { /* Assume an edge, not corner */ 64147c6ae99SBarry Smith n6 = rank - m*(n-1) + m-1 - m*n; 64247c6ae99SBarry Smith n15 = rank - m*(n-1) + m-1; 64347c6ae99SBarry Smith n24 = rank - m*(n-1) + m-1 + m*n; 64447c6ae99SBarry Smith } 64547c6ae99SBarry Smith 646ce00eea3SSatish Balay if ((xe==M) && (zs==0)) { /* Assume an edge, not corner */ 64747c6ae99SBarry Smith n2 = size - (m*n-rank) - (m-1) - m; 64847c6ae99SBarry Smith n5 = size - (m*n-rank) - (m-1); 64947c6ae99SBarry Smith n8 = size - (m*n-rank) - (m-1) + m; 65047c6ae99SBarry Smith } 65147c6ae99SBarry Smith 652ce00eea3SSatish Balay if ((xe==M) && (ze==P)) { /* Assume an edge, not corner */ 65347c6ae99SBarry Smith n20 = m*n - (size - rank) - (m-1) - m; 65447c6ae99SBarry Smith n23 = m*n - (size - rank) - (m-1); 65547c6ae99SBarry Smith n26 = m*n - (size - rank) - (m-1) + m; 65647c6ae99SBarry Smith } 65747c6ae99SBarry Smith 658ce00eea3SSatish Balay if ((xe==M) && (ys==0)) { /* Assume an edge, not corner */ 65947c6ae99SBarry Smith n2 = rank + m*(n-1) - (m-1) - m*n; 66047c6ae99SBarry Smith n11 = rank + m*(n-1) - (m-1); 66147c6ae99SBarry Smith n20 = rank + m*(n-1) - (m-1) + m*n; 66247c6ae99SBarry Smith } 66347c6ae99SBarry Smith 664ce00eea3SSatish Balay if ((xe==M) && (ye==N)) { /* Assume an edge, not corner */ 66547c6ae99SBarry Smith n8 = rank - m*n +1 - m*n; 66647c6ae99SBarry Smith n17 = rank - m*n +1; 66747c6ae99SBarry Smith n26 = rank - m*n +1 + m*n; 66847c6ae99SBarry Smith } 66947c6ae99SBarry Smith 67047c6ae99SBarry Smith if ((ys==0) && (zs==0)) { /* Assume an edge, not corner */ 67147c6ae99SBarry Smith n0 = size - m + rank -1; 67247c6ae99SBarry Smith n1 = size - m + rank; 67347c6ae99SBarry Smith n2 = size - m + rank +1; 67447c6ae99SBarry Smith } 67547c6ae99SBarry Smith 67647c6ae99SBarry Smith if ((ys==0) && (ze==P)) { /* Assume an edge, not corner */ 67747c6ae99SBarry Smith n18 = m*n - (size - rank) + m*(n-1) -1; 67847c6ae99SBarry Smith n19 = m*n - (size - rank) + m*(n-1); 67947c6ae99SBarry Smith n20 = m*n - (size - rank) + m*(n-1) +1; 68047c6ae99SBarry Smith } 68147c6ae99SBarry Smith 68247c6ae99SBarry Smith if ((ye==N) && (zs==0)) { /* Assume an edge, not corner */ 68347c6ae99SBarry Smith n6 = size - (m*n-rank) - m * (n-1) -1; 68447c6ae99SBarry Smith n7 = size - (m*n-rank) - m * (n-1); 68547c6ae99SBarry Smith n8 = size - (m*n-rank) - m * (n-1) +1; 68647c6ae99SBarry Smith } 68747c6ae99SBarry Smith 68847c6ae99SBarry Smith if ((ye==N) && (ze==P)) { /* Assume an edge, not corner */ 68947c6ae99SBarry Smith n24 = rank - (size-m) -1; 69047c6ae99SBarry Smith n25 = rank - (size-m); 69147c6ae99SBarry Smith n26 = rank - (size-m) +1; 69247c6ae99SBarry Smith } 69347c6ae99SBarry Smith 69447c6ae99SBarry Smith /* Check for Corners */ 6958865f1eaSKarl Rupp if ((xs==0) && (ys==0) && (zs==0)) n0 = size -1; 6968865f1eaSKarl Rupp if ((xs==0) && (ys==0) && (ze==P)) n18 = m*n-1; 6978865f1eaSKarl Rupp if ((xs==0) && (ye==N) && (zs==0)) n6 = (size-1)-m*(n-1); 6988865f1eaSKarl Rupp if ((xs==0) && (ye==N) && (ze==P)) n24 = m-1; 6998865f1eaSKarl Rupp if ((xe==M) && (ys==0) && (zs==0)) n2 = size-m; 7008865f1eaSKarl Rupp if ((xe==M) && (ys==0) && (ze==P)) n20 = m*n-m; 7018865f1eaSKarl Rupp if ((xe==M) && (ye==N) && (zs==0)) n8 = size-m*n; 7028865f1eaSKarl Rupp if ((xe==M) && (ye==N) && (ze==P)) n26 = 0; 70347c6ae99SBarry Smith 70447c6ae99SBarry Smith /* Check for when not X,Y, and Z Periodic */ 70547c6ae99SBarry Smith 70647c6ae99SBarry Smith /* If not X periodic */ 707bff4a2f0SMatthew G. Knepley if (bx != DM_BOUNDARY_PERIODIC) { 7088865f1eaSKarl Rupp if (xs==0) n0 = n3 = n6 = n9 = n12 = n15 = n18 = n21 = n24 = -2; 7098865f1eaSKarl Rupp if (xe==M) n2 = n5 = n8 = n11 = n14 = n17 = n20 = n23 = n26 = -2; 71047c6ae99SBarry Smith } 71147c6ae99SBarry Smith 71247c6ae99SBarry Smith /* If not Y periodic */ 713bff4a2f0SMatthew G. Knepley if (by != DM_BOUNDARY_PERIODIC) { 7148865f1eaSKarl Rupp if (ys==0) n0 = n1 = n2 = n9 = n10 = n11 = n18 = n19 = n20 = -2; 7158865f1eaSKarl Rupp if (ye==N) n6 = n7 = n8 = n15 = n16 = n17 = n24 = n25 = n26 = -2; 71647c6ae99SBarry Smith } 71747c6ae99SBarry Smith 71847c6ae99SBarry Smith /* If not Z periodic */ 719bff4a2f0SMatthew G. Knepley if (bz != DM_BOUNDARY_PERIODIC) { 7208865f1eaSKarl Rupp if (zs==0) n0 = n1 = n2 = n3 = n4 = n5 = n6 = n7 = n8 = -2; 7218865f1eaSKarl Rupp if (ze==P) n18 = n19 = n20 = n21 = n22 = n23 = n24 = n25 = n26 = -2; 72247c6ae99SBarry Smith } 72347c6ae99SBarry Smith 7249566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(27,&dd->neighbors)); 7258865f1eaSKarl Rupp 72647c6ae99SBarry Smith dd->neighbors[0] = n0; 72747c6ae99SBarry Smith dd->neighbors[1] = n1; 72847c6ae99SBarry Smith dd->neighbors[2] = n2; 72947c6ae99SBarry Smith dd->neighbors[3] = n3; 73047c6ae99SBarry Smith dd->neighbors[4] = n4; 73147c6ae99SBarry Smith dd->neighbors[5] = n5; 73247c6ae99SBarry Smith dd->neighbors[6] = n6; 73347c6ae99SBarry Smith dd->neighbors[7] = n7; 73447c6ae99SBarry Smith dd->neighbors[8] = n8; 73547c6ae99SBarry Smith dd->neighbors[9] = n9; 73647c6ae99SBarry Smith dd->neighbors[10] = n10; 73747c6ae99SBarry Smith dd->neighbors[11] = n11; 73847c6ae99SBarry Smith dd->neighbors[12] = n12; 73947c6ae99SBarry Smith dd->neighbors[13] = rank; 74047c6ae99SBarry Smith dd->neighbors[14] = n14; 74147c6ae99SBarry Smith dd->neighbors[15] = n15; 74247c6ae99SBarry Smith dd->neighbors[16] = n16; 74347c6ae99SBarry Smith dd->neighbors[17] = n17; 74447c6ae99SBarry Smith dd->neighbors[18] = n18; 74547c6ae99SBarry Smith dd->neighbors[19] = n19; 74647c6ae99SBarry Smith dd->neighbors[20] = n20; 74747c6ae99SBarry Smith dd->neighbors[21] = n21; 74847c6ae99SBarry Smith dd->neighbors[22] = n22; 74947c6ae99SBarry Smith dd->neighbors[23] = n23; 75047c6ae99SBarry Smith dd->neighbors[24] = n24; 75147c6ae99SBarry Smith dd->neighbors[25] = n25; 75247c6ae99SBarry Smith dd->neighbors[26] = n26; 75347c6ae99SBarry Smith 75447c6ae99SBarry Smith /* If star stencil then delete the corner neighbors */ 755d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_STAR) { 75647c6ae99SBarry Smith /* save information about corner neighbors */ 75747c6ae99SBarry Smith sn0 = n0; sn1 = n1; sn2 = n2; sn3 = n3; sn5 = n5; sn6 = n6; sn7 = n7; 75847c6ae99SBarry Smith sn8 = n8; sn9 = n9; sn11 = n11; sn15 = n15; sn17 = n17; sn18 = n18; 75947c6ae99SBarry Smith sn19 = n19; sn20 = n20; sn21 = n21; sn23 = n23; sn24 = n24; sn25 = n25; 76047c6ae99SBarry Smith sn26 = n26; 7618865f1eaSKarl Rupp n0 = n1 = n2 = n3 = n5 = n6 = n7 = n8 = n9 = n11 = n15 = n17 = n18 = n19 = n20 = n21 = n23 = n24 = n25 = n26 = -1; 76247c6ae99SBarry Smith } 76347c6ae99SBarry Smith 7649566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((Xe-Xs)*(Ye-Ys)*(Ze-Zs),&idx)); 76547c6ae99SBarry Smith 76647c6ae99SBarry Smith nn = 0; 76747c6ae99SBarry Smith /* Bottom Level */ 76847c6ae99SBarry Smith for (k=0; k<s_z; k++) { 76947c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 77047c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 771ce00eea3SSatish Balay x_t = lx[n0 % m]; 77247c6ae99SBarry Smith y_t = ly[(n0 % (m*n))/m]; 77347c6ae99SBarry Smith z_t = lz[n0 / (m*n)]; 77447c6ae99SBarry 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; 7758865f1eaSKarl 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 */ 7768865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 77747c6ae99SBarry Smith } 77847c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 77947c6ae99SBarry Smith x_t = x; 78047c6ae99SBarry Smith y_t = ly[(n1 % (m*n))/m]; 78147c6ae99SBarry Smith z_t = lz[n1 / (m*n)]; 78247c6ae99SBarry 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; 7838865f1eaSKarl Rupp if (twod && (s_t < 0)) s_t = bases[n1] + x_t*y_t*z_t - (s_y+1-i)*x_t; /* 2D case */ 7848865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 78547c6ae99SBarry Smith } 78647c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 787ce00eea3SSatish Balay x_t = lx[n2 % m]; 78847c6ae99SBarry Smith y_t = ly[(n2 % (m*n))/m]; 78947c6ae99SBarry Smith z_t = lz[n2 / (m*n)]; 79047c6ae99SBarry 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; 7918865f1eaSKarl Rupp if (twod && (s_t < 0)) s_t = bases[n2] + x_t*y_t*z_t - (s_y+1-i)*x_t; /* 2D case */ 7928865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 79347c6ae99SBarry Smith } 79447c6ae99SBarry Smith } 79547c6ae99SBarry Smith 79647c6ae99SBarry Smith for (i=0; i<y; i++) { 79747c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 798ce00eea3SSatish Balay x_t = lx[n3 % m]; 79947c6ae99SBarry Smith y_t = y; 80047c6ae99SBarry Smith z_t = lz[n3 / (m*n)]; 80147c6ae99SBarry Smith s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8028865f1eaSKarl 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 */ 8038865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 80447c6ae99SBarry Smith } 80547c6ae99SBarry Smith 80647c6ae99SBarry Smith if (n4 >= 0) { /* middle */ 80747c6ae99SBarry Smith x_t = x; 80847c6ae99SBarry Smith y_t = y; 80947c6ae99SBarry Smith z_t = lz[n4 / (m*n)]; 81047c6ae99SBarry Smith s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8118865f1eaSKarl 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 */ 8128865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 813bff4a2f0SMatthew G. Knepley } else if (bz == DM_BOUNDARY_MIRROR) { 814a6fd6b0aSBarry Smith for (j=0; j<x; j++) idx[nn++] = bases[rank] + j + x*i + (s_z - k - 1)*x*y; 81547c6ae99SBarry Smith } 81647c6ae99SBarry Smith 81747c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 818ce00eea3SSatish Balay x_t = lx[n5 % m]; 81947c6ae99SBarry Smith y_t = y; 82047c6ae99SBarry Smith z_t = lz[n5 / (m*n)]; 82147c6ae99SBarry Smith s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8228865f1eaSKarl 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 */ 8238865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 82447c6ae99SBarry Smith } 82547c6ae99SBarry Smith } 82647c6ae99SBarry Smith 82747c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 82847c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 829ce00eea3SSatish Balay x_t = lx[n6 % m]; 83047c6ae99SBarry Smith y_t = ly[(n6 % (m*n))/m]; 83147c6ae99SBarry Smith z_t = lz[n6 / (m*n)]; 83247c6ae99SBarry Smith s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8338865f1eaSKarl 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 */ 8348865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 83547c6ae99SBarry Smith } 83647c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 83747c6ae99SBarry Smith x_t = x; 83847c6ae99SBarry Smith y_t = ly[(n7 % (m*n))/m]; 83947c6ae99SBarry Smith z_t = lz[n7 / (m*n)]; 84047c6ae99SBarry Smith s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8418865f1eaSKarl 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 */ 8428865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 84347c6ae99SBarry Smith } 84447c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 845ce00eea3SSatish Balay x_t = lx[n8 % m]; 84647c6ae99SBarry Smith y_t = ly[(n8 % (m*n))/m]; 84747c6ae99SBarry Smith z_t = lz[n8 / (m*n)]; 84847c6ae99SBarry Smith s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8498865f1eaSKarl 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 */ 8508865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 85147c6ae99SBarry Smith } 85247c6ae99SBarry Smith } 85347c6ae99SBarry Smith } 85447c6ae99SBarry Smith 85547c6ae99SBarry Smith /* Middle Level */ 85647c6ae99SBarry Smith for (k=0; k<z; k++) { 85747c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 85847c6ae99SBarry Smith if (n9 >= 0) { /* left below */ 859ce00eea3SSatish Balay x_t = lx[n9 % m]; 86047c6ae99SBarry Smith y_t = ly[(n9 % (m*n))/m]; 86147c6ae99SBarry Smith /* z_t = z; */ 86247c6ae99SBarry Smith s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 8638865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 86447c6ae99SBarry Smith } 86547c6ae99SBarry Smith if (n10 >= 0) { /* directly below */ 86647c6ae99SBarry Smith x_t = x; 86747c6ae99SBarry Smith y_t = ly[(n10 % (m*n))/m]; 86847c6ae99SBarry Smith /* z_t = z; */ 86947c6ae99SBarry Smith s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 8708865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 871bff4a2f0SMatthew G. Knepley } else if (by == DM_BOUNDARY_MIRROR) { 872a6fd6b0aSBarry Smith for (j=0; j<x; j++) idx[nn++] = bases[rank] + j + k*x*y + (s_y - i)*x; 87347c6ae99SBarry Smith } 87447c6ae99SBarry Smith if (n11 >= 0) { /* right below */ 875ce00eea3SSatish Balay x_t = lx[n11 % m]; 87647c6ae99SBarry Smith y_t = ly[(n11 % (m*n))/m]; 87747c6ae99SBarry Smith /* z_t = z; */ 87847c6ae99SBarry Smith s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 8798865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 88047c6ae99SBarry Smith } 88147c6ae99SBarry Smith } 88247c6ae99SBarry Smith 88347c6ae99SBarry Smith for (i=0; i<y; i++) { 88447c6ae99SBarry Smith if (n12 >= 0) { /* directly left */ 885ce00eea3SSatish Balay x_t = lx[n12 % m]; 88647c6ae99SBarry Smith y_t = y; 88747c6ae99SBarry Smith /* z_t = z; */ 88847c6ae99SBarry Smith s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t; 8898865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 890bff4a2f0SMatthew G. Knepley } else if (bx == DM_BOUNDARY_MIRROR) { 891a6fd6b0aSBarry Smith for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + s_x - j - 1 + k*x*y + i*x; 89247c6ae99SBarry Smith } 89347c6ae99SBarry Smith 89447c6ae99SBarry Smith /* Interior */ 89547c6ae99SBarry Smith s_t = bases[rank] + i*x + k*x*y; 8968865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = s_t++; 89747c6ae99SBarry Smith 89847c6ae99SBarry Smith if (n14 >= 0) { /* directly right */ 899ce00eea3SSatish Balay x_t = lx[n14 % m]; 90047c6ae99SBarry Smith y_t = y; 90147c6ae99SBarry Smith /* z_t = z; */ 90247c6ae99SBarry Smith s_t = bases[n14] + i*x_t + k*x_t*y_t; 9038865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 904bff4a2f0SMatthew G. Knepley } else if (bx == DM_BOUNDARY_MIRROR) { 905a6fd6b0aSBarry Smith for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x - j - 1 + k*x*y + i*x; 90647c6ae99SBarry Smith } 90747c6ae99SBarry Smith } 90847c6ae99SBarry Smith 90947c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 91047c6ae99SBarry Smith if (n15 >= 0) { /* left above */ 911ce00eea3SSatish Balay x_t = lx[n15 % m]; 91247c6ae99SBarry Smith y_t = ly[(n15 % (m*n))/m]; 91347c6ae99SBarry Smith /* z_t = z; */ 91447c6ae99SBarry Smith s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t; 9158865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 91647c6ae99SBarry Smith } 91747c6ae99SBarry Smith if (n16 >= 0) { /* directly above */ 91847c6ae99SBarry Smith x_t = x; 91947c6ae99SBarry Smith y_t = ly[(n16 % (m*n))/m]; 92047c6ae99SBarry Smith /* z_t = z; */ 92147c6ae99SBarry Smith s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t; 9228865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 923bff4a2f0SMatthew G. Knepley } else if (by == DM_BOUNDARY_MIRROR) { 924a6fd6b0aSBarry Smith for (j=0; j<x; j++) idx[nn++] = bases[rank] + j + k*x*y + (y-i)*x; 92547c6ae99SBarry Smith } 92647c6ae99SBarry Smith if (n17 >= 0) { /* right above */ 927ce00eea3SSatish Balay x_t = lx[n17 % m]; 92847c6ae99SBarry Smith y_t = ly[(n17 % (m*n))/m]; 92947c6ae99SBarry Smith /* z_t = z; */ 93047c6ae99SBarry Smith s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t; 9318865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 93247c6ae99SBarry Smith } 93347c6ae99SBarry Smith } 93447c6ae99SBarry Smith } 93547c6ae99SBarry Smith 93647c6ae99SBarry Smith /* Upper Level */ 93747c6ae99SBarry Smith for (k=0; k<s_z; k++) { 93847c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 93947c6ae99SBarry Smith if (n18 >= 0) { /* left below */ 940ce00eea3SSatish Balay x_t = lx[n18 % m]; 94147c6ae99SBarry Smith y_t = ly[(n18 % (m*n))/m]; 94247c6ae99SBarry Smith /* z_t = lz[n18 / (m*n)]; */ 94347c6ae99SBarry Smith s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 9448865f1eaSKarl 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 */ 9458865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 94647c6ae99SBarry Smith } 94747c6ae99SBarry Smith if (n19 >= 0) { /* directly below */ 94847c6ae99SBarry Smith x_t = x; 94947c6ae99SBarry Smith y_t = ly[(n19 % (m*n))/m]; 95047c6ae99SBarry Smith /* z_t = lz[n19 / (m*n)]; */ 95147c6ae99SBarry Smith s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 9528865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n19] - (s_y+1-i)*x_t + x_t*y_t; /* 2d case */ 9538865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 95447c6ae99SBarry Smith } 95547c6ae99SBarry Smith if (n20 >= 0) { /* right below */ 956ce00eea3SSatish Balay x_t = lx[n20 % m]; 95747c6ae99SBarry Smith y_t = ly[(n20 % (m*n))/m]; 95847c6ae99SBarry Smith /* z_t = lz[n20 / (m*n)]; */ 95947c6ae99SBarry Smith s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 9608865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n20] - (s_y+1-i)*x_t + x_t*y_t; /* 2d case */ 9618865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 96247c6ae99SBarry Smith } 96347c6ae99SBarry Smith } 96447c6ae99SBarry Smith 96547c6ae99SBarry Smith for (i=0; i<y; i++) { 96647c6ae99SBarry Smith if (n21 >= 0) { /* directly left */ 967ce00eea3SSatish Balay x_t = lx[n21 % m]; 96847c6ae99SBarry Smith y_t = y; 96947c6ae99SBarry Smith /* z_t = lz[n21 / (m*n)]; */ 97047c6ae99SBarry Smith s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t; 9718865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n21] + (i+1)*x_t - s_x; /* 2d case */ 9728865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 97347c6ae99SBarry Smith } 97447c6ae99SBarry Smith 97547c6ae99SBarry Smith if (n22 >= 0) { /* middle */ 97647c6ae99SBarry Smith x_t = x; 97747c6ae99SBarry Smith y_t = y; 97847c6ae99SBarry Smith /* z_t = lz[n22 / (m*n)]; */ 97947c6ae99SBarry Smith s_t = bases[n22] + i*x_t + k*x_t*y_t; 9808865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n22] + i*x_t; /* 2d case */ 9818865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 982bff4a2f0SMatthew G. Knepley } else if (bz == DM_BOUNDARY_MIRROR) { 983a6fd6b0aSBarry Smith for (j=0; j<x; j++) idx[nn++] = bases[rank] + j + (z-k-1)*x*y + i*x; 98447c6ae99SBarry Smith } 98547c6ae99SBarry Smith 98647c6ae99SBarry Smith if (n23 >= 0) { /* directly right */ 987ce00eea3SSatish Balay x_t = lx[n23 % m]; 98847c6ae99SBarry Smith y_t = y; 98947c6ae99SBarry Smith /* z_t = lz[n23 / (m*n)]; */ 99047c6ae99SBarry Smith s_t = bases[n23] + i*x_t + k*x_t*y_t; 9918865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n23] + i*x_t; /* 2d case */ 9928865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 99347c6ae99SBarry Smith } 99447c6ae99SBarry Smith } 99547c6ae99SBarry Smith 99647c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 99747c6ae99SBarry Smith if (n24 >= 0) { /* left above */ 998ce00eea3SSatish Balay x_t = lx[n24 % m]; 99947c6ae99SBarry Smith y_t = ly[(n24 % (m*n))/m]; 100047c6ae99SBarry Smith /* z_t = lz[n24 / (m*n)]; */ 100147c6ae99SBarry Smith s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t; 10028865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n24] + i*x_t - s_x; /* 2d case */ 10038865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 100447c6ae99SBarry Smith } 100547c6ae99SBarry Smith if (n25 >= 0) { /* directly above */ 100647c6ae99SBarry Smith x_t = x; 100747c6ae99SBarry Smith y_t = ly[(n25 % (m*n))/m]; 100847c6ae99SBarry Smith /* z_t = lz[n25 / (m*n)]; */ 100947c6ae99SBarry Smith s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t; 10108865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n25] + (i-1)*x_t; /* 2d case */ 10118865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 101247c6ae99SBarry Smith } 101347c6ae99SBarry Smith if (n26 >= 0) { /* right above */ 1014ce00eea3SSatish Balay x_t = lx[n26 % m]; 101547c6ae99SBarry Smith y_t = ly[(n26 % (m*n))/m]; 101647c6ae99SBarry Smith /* z_t = lz[n26 / (m*n)]; */ 101747c6ae99SBarry Smith s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t; 10188865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n26] + (i-1)*x_t; /* 2d case */ 10198865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 102047c6ae99SBarry Smith } 102147c6ae99SBarry Smith } 102247c6ae99SBarry Smith } 1023ce00eea3SSatish Balay 10249566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(comm,dof,nn,idx,PETSC_USE_POINTER,&from)); 10259566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(global,from,local,to,>ol)); 10269566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)da,(PetscObject)gtol)); 10279566063dSJacob Faibussowitsch PetscCall(ISDestroy(&to)); 10289566063dSJacob Faibussowitsch PetscCall(ISDestroy(&from)); 102947c6ae99SBarry Smith 1030d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_STAR) { 103147c6ae99SBarry Smith n0 = sn0; n1 = sn1; n2 = sn2; n3 = sn3; n5 = sn5; n6 = sn6; n7 = sn7; 103247c6ae99SBarry Smith n8 = sn8; n9 = sn9; n11 = sn11; n15 = sn15; n17 = sn17; n18 = sn18; 103347c6ae99SBarry Smith n19 = sn19; n20 = sn20; n21 = sn21; n23 = sn23; n24 = sn24; n25 = sn25; 103447c6ae99SBarry Smith n26 = sn26; 1035ce00eea3SSatish Balay } 103647c6ae99SBarry Smith 1037288e7d53SBarry Smith if (((stencil_type == DMDA_STENCIL_STAR) || (bx != DM_BOUNDARY_PERIODIC && bx) || (by != DM_BOUNDARY_PERIODIC && by) || (bz != DM_BOUNDARY_PERIODIC && bz))) { 1038ce00eea3SSatish Balay /* 1039ce00eea3SSatish Balay Recompute the local to global mappings, this time keeping the 1040ce00eea3SSatish Balay information about the cross corner processor numbers. 1041ce00eea3SSatish Balay */ 104247c6ae99SBarry Smith nn = 0; 104347c6ae99SBarry Smith /* Bottom Level */ 104447c6ae99SBarry Smith for (k=0; k<s_z; k++) { 104547c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 104647c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 1047ce00eea3SSatish Balay x_t = lx[n0 % m]; 104847c6ae99SBarry Smith y_t = ly[(n0 % (m*n))/m]; 104947c6ae99SBarry Smith z_t = lz[n0 / (m*n)]; 105047c6ae99SBarry 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; 10518865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1052ce00eea3SSatish Balay } else if (Xs-xs < 0 && Ys-ys < 0 && Zs-zs < 0) { 10538865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 105447c6ae99SBarry Smith } 105547c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 105647c6ae99SBarry Smith x_t = x; 105747c6ae99SBarry Smith y_t = ly[(n1 % (m*n))/m]; 105847c6ae99SBarry Smith z_t = lz[n1 / (m*n)]; 105947c6ae99SBarry 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; 10608865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1061ce00eea3SSatish Balay } else if (Ys-ys < 0 && Zs-zs < 0) { 10628865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 106347c6ae99SBarry Smith } 106447c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 1065ce00eea3SSatish Balay x_t = lx[n2 % m]; 106647c6ae99SBarry Smith y_t = ly[(n2 % (m*n))/m]; 106747c6ae99SBarry Smith z_t = lz[n2 / (m*n)]; 106847c6ae99SBarry 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; 10698865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1070ce00eea3SSatish Balay } else if (xe-Xe < 0 && Ys-ys < 0 && Zs-zs < 0) { 10718865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 107247c6ae99SBarry Smith } 107347c6ae99SBarry Smith } 107447c6ae99SBarry Smith 107547c6ae99SBarry Smith for (i=0; i<y; i++) { 107647c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 1077ce00eea3SSatish Balay x_t = lx[n3 % m]; 107847c6ae99SBarry Smith y_t = y; 107947c6ae99SBarry Smith z_t = lz[n3 / (m*n)]; 108047c6ae99SBarry Smith s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 10818865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1082ce00eea3SSatish Balay } else if (Xs-xs < 0 && Zs-zs < 0) { 10838865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 108447c6ae99SBarry Smith } 108547c6ae99SBarry Smith 108647c6ae99SBarry Smith if (n4 >= 0) { /* middle */ 108747c6ae99SBarry Smith x_t = x; 108847c6ae99SBarry Smith y_t = y; 108947c6ae99SBarry Smith z_t = lz[n4 / (m*n)]; 109047c6ae99SBarry Smith s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 10918865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1092ce00eea3SSatish Balay } else if (Zs-zs < 0) { 1093bff4a2f0SMatthew G. Knepley if (bz == DM_BOUNDARY_MIRROR) { 1094a6fd6b0aSBarry Smith for (j=0; j<x; j++) idx[nn++] = bases[rank] + j + x*i + (s_z - k - 1)*x*y; 1095c2859e5eSBarry Smith } else { 10968865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 109747c6ae99SBarry Smith } 1098c2859e5eSBarry Smith } 109947c6ae99SBarry Smith 110047c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 1101ce00eea3SSatish Balay x_t = lx[n5 % m]; 110247c6ae99SBarry Smith y_t = y; 110347c6ae99SBarry Smith z_t = lz[n5 / (m*n)]; 110447c6ae99SBarry Smith s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 11058865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1106ce00eea3SSatish Balay } else if (xe-Xe < 0 && Zs-zs < 0) { 11078865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 110847c6ae99SBarry Smith } 110947c6ae99SBarry Smith } 111047c6ae99SBarry Smith 111147c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 111247c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 1113ce00eea3SSatish Balay x_t = lx[n6 % m]; 111447c6ae99SBarry Smith y_t = ly[(n6 % (m*n))/m]; 111547c6ae99SBarry Smith z_t = lz[n6 / (m*n)]; 111647c6ae99SBarry Smith s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 11178865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1118ce00eea3SSatish Balay } else if (Xs-xs < 0 && ye-Ye < 0 && Zs-zs < 0) { 11198865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 112047c6ae99SBarry Smith } 112147c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 112247c6ae99SBarry Smith x_t = x; 112347c6ae99SBarry Smith y_t = ly[(n7 % (m*n))/m]; 112447c6ae99SBarry Smith z_t = lz[n7 / (m*n)]; 112547c6ae99SBarry Smith s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 11268865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1127ce00eea3SSatish Balay } else if (ye-Ye < 0 && Zs-zs < 0) { 11288865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 112947c6ae99SBarry Smith } 113047c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 1131ce00eea3SSatish Balay x_t = lx[n8 % m]; 113247c6ae99SBarry Smith y_t = ly[(n8 % (m*n))/m]; 113347c6ae99SBarry Smith z_t = lz[n8 / (m*n)]; 113447c6ae99SBarry Smith s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 11358865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1136ce00eea3SSatish Balay } else if (xe-Xe < 0 && ye-Ye < 0 && Zs-zs < 0) { 11378865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 113847c6ae99SBarry Smith } 113947c6ae99SBarry Smith } 114047c6ae99SBarry Smith } 114147c6ae99SBarry Smith 114247c6ae99SBarry Smith /* Middle Level */ 114347c6ae99SBarry Smith for (k=0; k<z; k++) { 114447c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 114547c6ae99SBarry Smith if (n9 >= 0) { /* left below */ 1146ce00eea3SSatish Balay x_t = lx[n9 % m]; 114747c6ae99SBarry Smith y_t = ly[(n9 % (m*n))/m]; 114847c6ae99SBarry Smith /* z_t = z; */ 114947c6ae99SBarry Smith s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 11508865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1151ce00eea3SSatish Balay } else if (Xs-xs < 0 && Ys-ys < 0) { 11528865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 115347c6ae99SBarry Smith } 115447c6ae99SBarry Smith if (n10 >= 0) { /* directly below */ 115547c6ae99SBarry Smith x_t = x; 115647c6ae99SBarry Smith y_t = ly[(n10 % (m*n))/m]; 115747c6ae99SBarry Smith /* z_t = z; */ 115847c6ae99SBarry Smith s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 11598865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1160ce00eea3SSatish Balay } else if (Ys-ys < 0) { 1161bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_MIRROR) { 1162feb237baSPierre Jolivet for (j=0; j<x; j++) idx[nn++] = bases[rank] + j + k*x*y + (s_y - i)*x; 1163c2859e5eSBarry Smith } else { 11648865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 1165c2859e5eSBarry Smith } 116647c6ae99SBarry Smith } 116747c6ae99SBarry Smith if (n11 >= 0) { /* right below */ 1168ce00eea3SSatish Balay x_t = lx[n11 % m]; 116947c6ae99SBarry Smith y_t = ly[(n11 % (m*n))/m]; 117047c6ae99SBarry Smith /* z_t = z; */ 117147c6ae99SBarry Smith s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 11728865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1173ce00eea3SSatish Balay } else if (xe-Xe < 0 && Ys-ys < 0) { 11748865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 117547c6ae99SBarry Smith } 117647c6ae99SBarry Smith } 117747c6ae99SBarry Smith 117847c6ae99SBarry Smith for (i=0; i<y; i++) { 117947c6ae99SBarry Smith if (n12 >= 0) { /* directly left */ 1180ce00eea3SSatish Balay x_t = lx[n12 % m]; 118147c6ae99SBarry Smith y_t = y; 118247c6ae99SBarry Smith /* z_t = z; */ 118347c6ae99SBarry Smith s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t; 11848865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1185ce00eea3SSatish Balay } else if (Xs-xs < 0) { 1186bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_MIRROR) { 1187a6fd6b0aSBarry Smith for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + s_x - j - 1 + k*x*y + i*x; 1188c2859e5eSBarry Smith } else { 11898865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 119047c6ae99SBarry Smith } 1191c2859e5eSBarry Smith } 119247c6ae99SBarry Smith 119347c6ae99SBarry Smith /* Interior */ 119447c6ae99SBarry Smith s_t = bases[rank] + i*x + k*x*y; 11958865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = s_t++; 119647c6ae99SBarry Smith 119747c6ae99SBarry Smith if (n14 >= 0) { /* directly right */ 1198ce00eea3SSatish Balay x_t = lx[n14 % m]; 119947c6ae99SBarry Smith y_t = y; 120047c6ae99SBarry Smith /* z_t = z; */ 120147c6ae99SBarry Smith s_t = bases[n14] + i*x_t + k*x_t*y_t; 12028865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1203ce00eea3SSatish Balay } else if (xe-Xe < 0) { 1204bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_MIRROR) { 1205a6fd6b0aSBarry Smith for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x - j - 1 + k*x*y + i*x; 1206c2859e5eSBarry Smith } else { 12078865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 120847c6ae99SBarry Smith } 120947c6ae99SBarry Smith } 1210c2859e5eSBarry Smith } 121147c6ae99SBarry Smith 121247c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 121347c6ae99SBarry Smith if (n15 >= 0) { /* left above */ 1214ce00eea3SSatish Balay x_t = lx[n15 % m]; 121547c6ae99SBarry Smith y_t = ly[(n15 % (m*n))/m]; 121647c6ae99SBarry Smith /* z_t = z; */ 121747c6ae99SBarry Smith s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t; 12188865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1219ce00eea3SSatish Balay } else if (Xs-xs < 0 && ye-Ye < 0) { 12208865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 122147c6ae99SBarry Smith } 122247c6ae99SBarry Smith if (n16 >= 0) { /* directly above */ 122347c6ae99SBarry Smith x_t = x; 122447c6ae99SBarry Smith y_t = ly[(n16 % (m*n))/m]; 122547c6ae99SBarry Smith /* z_t = z; */ 122647c6ae99SBarry Smith s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t; 12278865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1228ce00eea3SSatish Balay } else if (ye-Ye < 0) { 1229bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_MIRROR) { 1230a6fd6b0aSBarry Smith for (j=0; j<x; j++) idx[nn++] = bases[rank] + j + k*x*y + (y-i)*x; 1231c2859e5eSBarry Smith } else { 12328865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 123347c6ae99SBarry Smith } 1234c2859e5eSBarry Smith } 123547c6ae99SBarry Smith if (n17 >= 0) { /* right above */ 1236ce00eea3SSatish Balay x_t = lx[n17 % m]; 123747c6ae99SBarry Smith y_t = ly[(n17 % (m*n))/m]; 123847c6ae99SBarry Smith /* z_t = z; */ 123947c6ae99SBarry Smith s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t; 12408865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1241ce00eea3SSatish Balay } else if (xe-Xe < 0 && ye-Ye < 0) { 12428865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 124347c6ae99SBarry Smith } 124447c6ae99SBarry Smith } 124547c6ae99SBarry Smith } 124647c6ae99SBarry Smith 124747c6ae99SBarry Smith /* Upper Level */ 124847c6ae99SBarry Smith for (k=0; k<s_z; k++) { 124947c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 125047c6ae99SBarry Smith if (n18 >= 0) { /* left below */ 1251ce00eea3SSatish Balay x_t = lx[n18 % m]; 125247c6ae99SBarry Smith y_t = ly[(n18 % (m*n))/m]; 125347c6ae99SBarry Smith /* z_t = lz[n18 / (m*n)]; */ 125447c6ae99SBarry Smith s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 12558865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1256ce00eea3SSatish Balay } else if (Xs-xs < 0 && Ys-ys < 0 && ze-Ze < 0) { 12578865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 125847c6ae99SBarry Smith } 125947c6ae99SBarry Smith if (n19 >= 0) { /* directly below */ 126047c6ae99SBarry Smith x_t = x; 126147c6ae99SBarry Smith y_t = ly[(n19 % (m*n))/m]; 126247c6ae99SBarry Smith /* z_t = lz[n19 / (m*n)]; */ 126347c6ae99SBarry Smith s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 12648865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1265ce00eea3SSatish Balay } else if (Ys-ys < 0 && ze-Ze < 0) { 12668865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 126747c6ae99SBarry Smith } 126847c6ae99SBarry Smith if (n20 >= 0) { /* right below */ 1269ce00eea3SSatish Balay x_t = lx[n20 % m]; 127047c6ae99SBarry Smith y_t = ly[(n20 % (m*n))/m]; 127147c6ae99SBarry Smith /* z_t = lz[n20 / (m*n)]; */ 127247c6ae99SBarry Smith s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 12738865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1274ce00eea3SSatish Balay } else if (xe-Xe < 0 && Ys-ys < 0 && ze-Ze < 0) { 12758865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 127647c6ae99SBarry Smith } 127747c6ae99SBarry Smith } 127847c6ae99SBarry Smith 127947c6ae99SBarry Smith for (i=0; i<y; i++) { 128047c6ae99SBarry Smith if (n21 >= 0) { /* directly left */ 1281ce00eea3SSatish Balay x_t = lx[n21 % m]; 128247c6ae99SBarry Smith y_t = y; 128347c6ae99SBarry Smith /* z_t = lz[n21 / (m*n)]; */ 128447c6ae99SBarry Smith s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t; 12858865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1286ce00eea3SSatish Balay } else if (Xs-xs < 0 && ze-Ze < 0) { 12878865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 128847c6ae99SBarry Smith } 128947c6ae99SBarry Smith 129047c6ae99SBarry Smith if (n22 >= 0) { /* middle */ 129147c6ae99SBarry Smith x_t = x; 129247c6ae99SBarry Smith y_t = y; 129347c6ae99SBarry Smith /* z_t = lz[n22 / (m*n)]; */ 129447c6ae99SBarry Smith s_t = bases[n22] + i*x_t + k*x_t*y_t; 12958865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1296ce00eea3SSatish Balay } else if (ze-Ze < 0) { 1297bff4a2f0SMatthew G. Knepley if (bz == DM_BOUNDARY_MIRROR) { 1298a6fd6b0aSBarry Smith for (j=0; j<x; j++) idx[nn++] = bases[rank] + j + (z-k-1)*x*y + i*x; 1299c2859e5eSBarry Smith } else { 13008865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 130147c6ae99SBarry Smith } 1302c2859e5eSBarry Smith } 130347c6ae99SBarry Smith 130447c6ae99SBarry Smith if (n23 >= 0) { /* directly right */ 1305ce00eea3SSatish Balay x_t = lx[n23 % m]; 130647c6ae99SBarry Smith y_t = y; 130747c6ae99SBarry Smith /* z_t = lz[n23 / (m*n)]; */ 130847c6ae99SBarry Smith s_t = bases[n23] + i*x_t + k*x_t*y_t; 13098865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1310ce00eea3SSatish Balay } else if (xe-Xe < 0 && ze-Ze < 0) { 13118865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 131247c6ae99SBarry Smith } 131347c6ae99SBarry Smith } 131447c6ae99SBarry Smith 131547c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 131647c6ae99SBarry Smith if (n24 >= 0) { /* left above */ 1317ce00eea3SSatish Balay x_t = lx[n24 % m]; 131847c6ae99SBarry Smith y_t = ly[(n24 % (m*n))/m]; 131947c6ae99SBarry Smith /* z_t = lz[n24 / (m*n)]; */ 132047c6ae99SBarry Smith s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t; 13218865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1322ce00eea3SSatish Balay } else if (Xs-xs < 0 && ye-Ye < 0 && ze-Ze < 0) { 13238865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 132447c6ae99SBarry Smith } 132547c6ae99SBarry Smith if (n25 >= 0) { /* directly above */ 132647c6ae99SBarry Smith x_t = x; 132747c6ae99SBarry Smith y_t = ly[(n25 % (m*n))/m]; 132847c6ae99SBarry Smith /* z_t = lz[n25 / (m*n)]; */ 132947c6ae99SBarry Smith s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t; 13308865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1331ce00eea3SSatish Balay } else if (ye-Ye < 0 && ze-Ze < 0) { 13328865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 133347c6ae99SBarry Smith } 133447c6ae99SBarry Smith if (n26 >= 0) { /* right above */ 1335ce00eea3SSatish Balay x_t = lx[n26 % m]; 133647c6ae99SBarry Smith y_t = ly[(n26 % (m*n))/m]; 133747c6ae99SBarry Smith /* z_t = lz[n26 / (m*n)]; */ 133847c6ae99SBarry Smith s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t; 13398865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1340ce00eea3SSatish Balay } else if (xe-Xe < 0 && ye-Ye < 0 && ze-Ze < 0) { 13418865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 134247c6ae99SBarry Smith } 134347c6ae99SBarry Smith } 134447c6ae99SBarry Smith } 134547c6ae99SBarry Smith } 134647c6ae99SBarry Smith /* 134747c6ae99SBarry Smith Set the local to global ordering in the global vector, this allows use 134847c6ae99SBarry Smith of VecSetValuesLocal(). 134947c6ae99SBarry Smith */ 13509566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(comm,dof,nn,idx,PETSC_OWN_POINTER,&da->ltogmap)); 13519566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)da,(PetscObject)da->ltogmap)); 135247c6ae99SBarry Smith 13539566063dSJacob Faibussowitsch PetscCall(PetscFree2(bases,ldims)); 1354ce00eea3SSatish Balay dd->m = m; dd->n = n; dd->p = p; 1355ce00eea3SSatish Balay /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */ 1356ce00eea3SSatish Balay dd->xs = xs*dof; dd->xe = xe*dof; dd->ys = ys; dd->ye = ye; dd->zs = zs; dd->ze = ze; 1357ce00eea3SSatish Balay dd->Xs = Xs*dof; dd->Xe = Xe*dof; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = Zs; dd->Ze = Ze; 1358ce00eea3SSatish Balay 13599566063dSJacob Faibussowitsch PetscCall(VecDestroy(&local)); 13609566063dSJacob Faibussowitsch PetscCall(VecDestroy(&global)); 1361ce00eea3SSatish Balay 1362ce00eea3SSatish Balay dd->gtol = gtol; 1363ce00eea3SSatish Balay dd->base = base; 1364ce00eea3SSatish Balay da->ops->view = DMView_DA_3d; 13650298fd71SBarry Smith dd->ltol = NULL; 13660298fd71SBarry Smith dd->ao = NULL; 136747c6ae99SBarry Smith PetscFunctionReturn(0); 136847c6ae99SBarry Smith } 1369564755cdSBarry Smith 137047c6ae99SBarry Smith /*@C 1371aa219208SBarry Smith DMDACreate3d - Creates an object that will manage the communication of three-dimensional 137247c6ae99SBarry Smith regular array data that is distributed across some processors. 137347c6ae99SBarry Smith 1374d083f849SBarry Smith Collective 137547c6ae99SBarry Smith 137647c6ae99SBarry Smith Input Parameters: 137747c6ae99SBarry Smith + comm - MPI communicator 13781321219cSEthan Coon . bx,by,bz - type of ghost nodes the array have. 1379bff4a2f0SMatthew G. Knepley Use one of DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_PERIODIC. 1380aa219208SBarry Smith . stencil_type - Type of stencil (DMDA_STENCIL_STAR or DMDA_STENCIL_BOX) 1381897f7067SBarry Smith . M,N,P - global dimension in each direction of the array 138247c6ae99SBarry Smith . m,n,p - corresponding number of processors in each dimension 138347c6ae99SBarry Smith (or PETSC_DECIDE to have calculated) 138447c6ae99SBarry Smith . dof - number of degrees of freedom per node 138510d7c030SMatthew G Knepley . s - stencil width 138610d7c030SMatthew G Knepley - lx, ly, lz - arrays containing the number of nodes in each cell along 13870298fd71SBarry Smith the x, y, and z coordinates, or NULL. If non-null, these 138847c6ae99SBarry Smith must be of length as m,n,p and the corresponding 138947c6ae99SBarry Smith m,n, or p cannot be PETSC_DECIDE. Sum of the lx[] entries must be M, sum of 139047c6ae99SBarry Smith the ly[] must N, sum of the lz[] must be P 139147c6ae99SBarry Smith 139247c6ae99SBarry Smith Output Parameter: 139347c6ae99SBarry Smith . da - the resulting distributed array object 139447c6ae99SBarry Smith 139547c6ae99SBarry Smith Options Database Key: 1396706a11cbSBarry Smith + -dm_view - Calls DMView() at the conclusion of DMDACreate3d() 1397e43dc8daSBarry Smith . -da_grid_x <nx> - number of grid points in x direction 1398e43dc8daSBarry Smith . -da_grid_y <ny> - number of grid points in y direction 1399e43dc8daSBarry Smith . -da_grid_z <nz> - number of grid points in z direction 140047c6ae99SBarry Smith . -da_processors_x <MX> - number of processors in x direction 140147c6ae99SBarry Smith . -da_processors_y <MY> - number of processors in y direction 140247c6ae99SBarry Smith . -da_processors_z <MZ> - number of processors in z direction 1403e0f5d30fSBarry Smith . -da_refine_x <rx> - refinement ratio in x direction 1404e0f5d30fSBarry Smith . -da_refine_y <ry> - refinement ratio in y direction 1405e0f5d30fSBarry Smith . -da_refine_z <rz>- refinement ratio in z directio 1406e43dc8daSBarry Smith - -da_refine <n> - refine the DMDA n times before creating it 140747c6ae99SBarry Smith 140847c6ae99SBarry Smith Level: beginner 140947c6ae99SBarry Smith 141047c6ae99SBarry Smith Notes: 1411aa219208SBarry Smith The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the 1412aa219208SBarry Smith standard 7-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes 141347c6ae99SBarry Smith the standard 27-pt stencil. 141447c6ae99SBarry Smith 1415aa219208SBarry Smith The array data itself is NOT stored in the DMDA, it is stored in Vec objects; 1416564755cdSBarry Smith The appropriate vector objects can be obtained with calls to DMCreateGlobalVector() 1417564755cdSBarry Smith and DMCreateLocalVector() and calls to VecDuplicate() if more are needed. 141847c6ae99SBarry Smith 1419897f7067SBarry Smith You must call DMSetUp() after this call before using this DM. 1420897f7067SBarry Smith 1421897f7067SBarry Smith If you wish to use the options database to change values in the DMDA call DMSetFromOptions() after this call 1422897f7067SBarry Smith but before DMSetUp(). 1423897f7067SBarry Smith 1424db781477SPatrick Sanan .seealso: `DMDestroy()`, `DMView()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMGlobalToLocalBegin()`, `DMDAGetRefinementFactor()`, 1425db781477SPatrick Sanan `DMGlobalToLocalEnd()`, `DMLocalToGlobalBegin()`, `DMLocalToLocalBegin()`, `DMLocalToLocalEnd()`, `DMDASetRefinementFactor()`, 1426db781477SPatrick Sanan `DMDAGetInfo()`, `DMCreateGlobalVector()`, `DMCreateLocalVector()`, `DMDACreateNaturalVector()`, `DMLoad()`, `DMDAGetOwnershipRanges()`, 1427db781477SPatrick Sanan `DMStagCreate3d()` 142847c6ae99SBarry Smith 142947c6ae99SBarry Smith @*/ 1430bff4a2f0SMatthew G. Knepley PetscErrorCode DMDACreate3d(MPI_Comm comm,DMBoundaryType bx,DMBoundaryType by,DMBoundaryType bz,DMDAStencilType stencil_type,PetscInt M, 14319a42bb27SBarry 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) 143247c6ae99SBarry Smith { 143347c6ae99SBarry Smith PetscFunctionBegin; 14349566063dSJacob Faibussowitsch PetscCall(DMDACreate(comm, da)); 14359566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*da, 3)); 14369566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(*da, M, N, P)); 14379566063dSJacob Faibussowitsch PetscCall(DMDASetNumProcs(*da, m, n, p)); 14389566063dSJacob Faibussowitsch PetscCall(DMDASetBoundaryType(*da, bx, by, bz)); 14399566063dSJacob Faibussowitsch PetscCall(DMDASetDof(*da, dof)); 14409566063dSJacob Faibussowitsch PetscCall(DMDASetStencilType(*da, stencil_type)); 14419566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(*da, s)); 14429566063dSJacob Faibussowitsch PetscCall(DMDASetOwnershipRanges(*da, lx, ly, lz)); 144347c6ae99SBarry Smith PetscFunctionReturn(0); 144447c6ae99SBarry Smith } 1445