xref: /petsc/src/dm/impls/da/da2.c (revision 76bd3646d776ac16fc835ea742fa097f15759704)
19a42bb27SBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h>    /*I   "petscdmda.h"   I*/
39804daf3SBarry Smith #include <petscdraw.h>
447c6ae99SBarry Smith 
5e0877f53SBarry Smith static PetscErrorCode DMView_DA_2d(DM da,PetscViewer viewer)
647c6ae99SBarry Smith {
747c6ae99SBarry Smith   PetscErrorCode ierr;
847c6ae99SBarry Smith   PetscMPIInt    rank;
9c9493c35SStefano Zampini   PetscBool      iascii,isdraw,isglvis,isbinary;
1047c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
119a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
129a42bb27SBarry Smith   PetscBool ismatlab;
139a42bb27SBarry Smith #endif
1447c6ae99SBarry Smith 
1547c6ae99SBarry Smith   PetscFunctionBegin;
16ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)da),&rank);CHKERRQ(ierr);
1747c6ae99SBarry Smith 
18251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
19251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
20c9493c35SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERGLVIS,&isglvis);CHKERRQ(ierr);
21251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
229a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
23251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERMATLAB,&ismatlab);CHKERRQ(ierr);
249a42bb27SBarry Smith #endif
2547c6ae99SBarry Smith   if (iascii) {
2647c6ae99SBarry Smith     PetscViewerFormat format;
2747c6ae99SBarry Smith 
2847c6ae99SBarry Smith     ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
2976a8abe0SBarry Smith     if (format == PETSC_VIEWER_LOAD_BALANCE) {
3076a8abe0SBarry Smith       PetscInt      i,nmax = 0,nmin = PETSC_MAX_INT,navg = 0,*nz,nzlocal;
3176a8abe0SBarry Smith       DMDALocalInfo info;
3276a8abe0SBarry Smith       PetscMPIInt   size;
3376a8abe0SBarry Smith       ierr = MPI_Comm_size(PetscObjectComm((PetscObject)da),&size);CHKERRQ(ierr);
3476a8abe0SBarry Smith       ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
3576a8abe0SBarry Smith       nzlocal = info.xm*info.ym;
3676a8abe0SBarry Smith       ierr = PetscMalloc1(size,&nz);CHKERRQ(ierr);
3776a8abe0SBarry Smith       ierr = MPI_Allgather(&nzlocal,1,MPIU_INT,nz,1,MPIU_INT,PetscObjectComm((PetscObject)da));CHKERRQ(ierr);
3876a8abe0SBarry Smith       for (i=0; i<(PetscInt)size; i++) {
3976a8abe0SBarry Smith         nmax = PetscMax(nmax,nz[i]);
4076a8abe0SBarry Smith         nmin = PetscMin(nmin,nz[i]);
4176a8abe0SBarry Smith         navg += nz[i];
4276a8abe0SBarry Smith       }
4376a8abe0SBarry Smith       ierr = PetscFree(nz);CHKERRQ(ierr);
4476a8abe0SBarry Smith       navg = navg/size;
4576a8abe0SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Load Balance - Grid Points: Min %D  avg %D  max %D\n",nmin,navg,nmax);CHKERRQ(ierr);
4676a8abe0SBarry Smith       PetscFunctionReturn(0);
4776a8abe0SBarry Smith     }
488135c375SStefano Zampini     if (format != PETSC_VIEWER_ASCII_VTK && format != PETSC_VIEWER_ASCII_VTK_CELL && format != PETSC_VIEWER_ASCII_GLVIS) {
49aa219208SBarry Smith       DMDALocalInfo info;
50aa219208SBarry Smith       ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
511575c14dSBarry Smith       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
5247c6ae99SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Processor [%d] M %D N %D m %D n %D w %D s %D\n",rank,dd->M,dd->N,dd->m,dd->n,dd->w,dd->s);CHKERRQ(ierr);
5347c6ae99SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"X range of indices: %D %D, Y range of indices: %D %D\n",info.xs,info.xs+info.xm,info.ys,info.ys+info.ym);CHKERRQ(ierr);
5447c6ae99SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
551575c14dSBarry Smith       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
568135c375SStefano Zampini     } else if (format == PETSC_VIEWER_ASCII_GLVIS) {
578135c375SStefano Zampini       ierr = DMView_DA_GLVis(da,viewer);CHKERRQ(ierr);
583da9ae13SJed Brown     } else {
593da9ae13SJed Brown       ierr = DMView_DA_VTK(da,viewer);CHKERRQ(ierr);
6047c6ae99SBarry Smith     }
6147c6ae99SBarry Smith   } else if (isdraw) {
6247c6ae99SBarry Smith     PetscDraw      draw;
6347c6ae99SBarry Smith     double         ymin = -1*dd->s-1,ymax = dd->N+dd->s;
6447c6ae99SBarry Smith     double         xmin = -1*dd->s-1,xmax = dd->M+dd->s;
6547c6ae99SBarry Smith     double         x,y;
668ea3bf28SBarry Smith     PetscInt       base;
678ea3bf28SBarry Smith     const PetscInt *idx;
6847c6ae99SBarry Smith     char           node[10];
6947c6ae99SBarry Smith     PetscBool      isnull;
7047c6ae99SBarry Smith 
7147c6ae99SBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
725b399a63SLisandro Dalcin     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
735b399a63SLisandro Dalcin     if (isnull) PetscFunctionReturn(0);
7447c6ae99SBarry Smith 
755b399a63SLisandro Dalcin     ierr = PetscDrawCheckResizedWindow(draw);CHKERRQ(ierr);
765b399a63SLisandro Dalcin     ierr = PetscDrawClear(draw);CHKERRQ(ierr);
775b399a63SLisandro Dalcin     ierr = PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);CHKERRQ(ierr);
785b399a63SLisandro Dalcin 
795b399a63SLisandro Dalcin     ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr);
8047c6ae99SBarry Smith     /* first processor draw all node lines */
8147c6ae99SBarry Smith     if (!rank) {
8247c6ae99SBarry Smith       ymin = 0.0; ymax = dd->N - 1;
8347c6ae99SBarry Smith       for (xmin=0; xmin<dd->M; xmin++) {
8447c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr);
8547c6ae99SBarry Smith       }
8647c6ae99SBarry Smith       xmin = 0.0; xmax = dd->M - 1;
8747c6ae99SBarry Smith       for (ymin=0; ymin<dd->N; ymin++) {
8847c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);CHKERRQ(ierr);
8947c6ae99SBarry Smith       }
9047c6ae99SBarry Smith     }
915b399a63SLisandro Dalcin     ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr);
925b399a63SLisandro Dalcin     ierr = PetscDrawFlush(draw);CHKERRQ(ierr);
9347c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
9447c6ae99SBarry Smith 
955b399a63SLisandro Dalcin     ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr);
9647c6ae99SBarry Smith     /* draw my box */
975b399a63SLisandro Dalcin     xmin = dd->xs/dd->w; xmax =(dd->xe-1)/dd->w; ymin = dd->ys; ymax = dd->ye - 1;
9847c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);CHKERRQ(ierr);
9947c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
10047c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
10147c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
10247c6ae99SBarry Smith     /* put in numbers */
10347c6ae99SBarry Smith     base = (dd->base)/dd->w;
10447c6ae99SBarry Smith     for (y=ymin; y<=ymax; y++) {
10547c6ae99SBarry Smith       for (x=xmin; x<=xmax; x++) {
1065b399a63SLisandro Dalcin         ierr = PetscSNPrintf(node,sizeof(node),"%d",(int)base++);CHKERRQ(ierr);
10747c6ae99SBarry Smith         ierr = PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node);CHKERRQ(ierr);
10847c6ae99SBarry Smith       }
10947c6ae99SBarry Smith     }
1105b399a63SLisandro Dalcin     ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr);
1115b399a63SLisandro Dalcin     ierr = PetscDrawFlush(draw);CHKERRQ(ierr);
11247c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
11347c6ae99SBarry Smith 
1145b399a63SLisandro Dalcin     ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr);
1155b399a63SLisandro Dalcin     /* overlay ghost numbers, useful for error checking */
11645b6f7e9SBarry Smith     ierr = ISLocalToGlobalMappingGetBlockIndices(da->ltogmap,&idx);CHKERRQ(ierr);
1175b399a63SLisandro Dalcin     base = 0; xmin = dd->Xs; xmax = dd->Xe; ymin = dd->Ys; ymax = dd->Ye;
11847c6ae99SBarry Smith     for (y=ymin; y<ymax; y++) {
11947c6ae99SBarry Smith       for (x=xmin; x<xmax; x++) {
12047c6ae99SBarry Smith         if ((base % dd->w) == 0) {
1215b399a63SLisandro Dalcin           ierr = PetscSNPrintf(node,sizeof(node),"%d",(int)(idx[base/dd->w]));CHKERRQ(ierr);
12247c6ae99SBarry Smith           ierr = PetscDrawString(draw,x/dd->w,y,PETSC_DRAW_BLUE,node);CHKERRQ(ierr);
12347c6ae99SBarry Smith         }
12447c6ae99SBarry Smith         base++;
12547c6ae99SBarry Smith       }
12647c6ae99SBarry Smith     }
127302440fdSBarry Smith     ierr = ISLocalToGlobalMappingRestoreBlockIndices(da->ltogmap,&idx);CHKERRQ(ierr);
1285b399a63SLisandro Dalcin     ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr);
1295b399a63SLisandro Dalcin     ierr = PetscDrawFlush(draw);CHKERRQ(ierr);
13047c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
131832b7cebSLisandro Dalcin     ierr = PetscDrawSave(draw);CHKERRQ(ierr);
132c9493c35SStefano Zampini   } else if (isglvis) {
133c9493c35SStefano Zampini     ierr = DMView_DA_GLVis(da,viewer);CHKERRQ(ierr);
1349a42bb27SBarry Smith   } else if (isbinary) {
1359a42bb27SBarry Smith     ierr = DMView_DA_Binary(da,viewer);CHKERRQ(ierr);
1369a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
1379a42bb27SBarry Smith   } else if (ismatlab) {
1389a42bb27SBarry Smith     ierr = DMView_DA_Matlab(da,viewer);CHKERRQ(ierr);
1399a42bb27SBarry Smith #endif
14011aeaf0aSBarry Smith   }
14147c6ae99SBarry Smith   PetscFunctionReturn(0);
14247c6ae99SBarry Smith }
14347c6ae99SBarry Smith 
14447c6ae99SBarry Smith 
14547c6ae99SBarry Smith #if defined(new)
14647c6ae99SBarry Smith /*
147aa219208SBarry Smith   DMDAGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix where local
148aa219208SBarry Smith     function lives on a DMDA
14947c6ae99SBarry Smith 
15047c6ae99SBarry Smith         y ~= (F(u + ha) - F(u))/h,
15147c6ae99SBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
15247c6ae99SBarry Smith         u = current iterate
15347c6ae99SBarry Smith         h = difference interval
15447c6ae99SBarry Smith */
155aa219208SBarry Smith PetscErrorCode DMDAGetDiagonal_MFFD(DM da,Vec U,Vec a)
15647c6ae99SBarry Smith {
15747c6ae99SBarry Smith   PetscScalar    h,*aa,*ww,v;
15847c6ae99SBarry Smith   PetscReal      epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
15947c6ae99SBarry Smith   PetscErrorCode ierr;
16047c6ae99SBarry Smith   PetscInt       gI,nI;
16147c6ae99SBarry Smith   MatStencil     stencil;
162aa219208SBarry Smith   DMDALocalInfo  info;
16347c6ae99SBarry Smith 
16447c6ae99SBarry Smith   PetscFunctionBegin;
16547c6ae99SBarry Smith   ierr = (*ctx->func)(0,U,a,ctx->funcctx);CHKERRQ(ierr);
16647c6ae99SBarry Smith   ierr = (*ctx->funcisetbase)(U,ctx->funcctx);CHKERRQ(ierr);
16747c6ae99SBarry Smith 
16847c6ae99SBarry Smith   ierr = VecGetArray(U,&ww);CHKERRQ(ierr);
16947c6ae99SBarry Smith   ierr = VecGetArray(a,&aa);CHKERRQ(ierr);
17047c6ae99SBarry Smith 
17147c6ae99SBarry Smith   nI = 0;
17247c6ae99SBarry Smith   h  = ww[gI];
17347c6ae99SBarry Smith   if (h == 0.0) h = 1.0;
17447c6ae99SBarry Smith   if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin;
17547c6ae99SBarry Smith   else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
17647c6ae99SBarry Smith   h *= epsilon;
17747c6ae99SBarry Smith 
17847c6ae99SBarry Smith   ww[gI] += h;
17947c6ae99SBarry Smith   ierr    = (*ctx->funci)(i,w,&v,ctx->funcctx);CHKERRQ(ierr);
18047c6ae99SBarry Smith   aa[nI]  = (v - aa[nI])/h;
18147c6ae99SBarry Smith   ww[gI] -= h;
18247c6ae99SBarry Smith   nI++;
1838865f1eaSKarl Rupp 
18447c6ae99SBarry Smith   ierr = VecRestoreArray(U,&ww);CHKERRQ(ierr);
18547c6ae99SBarry Smith   ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr);
18647c6ae99SBarry Smith   PetscFunctionReturn(0);
18747c6ae99SBarry Smith }
18847c6ae99SBarry Smith #endif
18947c6ae99SBarry Smith 
1907087cfbeSBarry Smith PetscErrorCode  DMSetUp_DA_2D(DM da)
19147c6ae99SBarry Smith {
19247c6ae99SBarry Smith   DM_DA            *dd = (DM_DA*)da->data;
19347c6ae99SBarry Smith   const PetscInt   M            = dd->M;
19447c6ae99SBarry Smith   const PetscInt   N            = dd->N;
19547c6ae99SBarry Smith   PetscInt         m            = dd->m;
19647c6ae99SBarry Smith   PetscInt         n            = dd->n;
19747c6ae99SBarry Smith   const PetscInt   dof          = dd->w;
19847c6ae99SBarry Smith   const PetscInt   s            = dd->s;
199bff4a2f0SMatthew G. Knepley   DMBoundaryType   bx           = dd->bx;
200bff4a2f0SMatthew G. Knepley   DMBoundaryType   by           = dd->by;
20119fd82e9SBarry Smith   DMDAStencilType  stencil_type = dd->stencil_type;
20247c6ae99SBarry Smith   PetscInt         *lx          = dd->lx;
20347c6ae99SBarry Smith   PetscInt         *ly          = dd->ly;
20447c6ae99SBarry Smith   MPI_Comm         comm;
20547c6ae99SBarry Smith   PetscMPIInt      rank,size;
206bd1fc5aeSBarry Smith   PetscInt         xs,xe,ys,ye,x,y,Xs,Xe,Ys,Ye,IXs,IXe,IYs,IYe;
2078ea3bf28SBarry Smith   PetscInt         up,down,left,right,i,n0,n1,n2,n3,n5,n6,n7,n8,*idx,nn;
208db87c5ecSEthan Coon   PetscInt         xbase,*bases,*ldims,j,x_t,y_t,s_t,base,count;
20947c6ae99SBarry Smith   PetscInt         s_x,s_y; /* s proportionalized to w */
21047c6ae99SBarry Smith   PetscInt         sn0 = 0,sn2 = 0,sn6 = 0,sn8 = 0;
21147c6ae99SBarry Smith   Vec              local,global;
212bd1fc5aeSBarry Smith   VecScatter       gtol;
21345b6f7e9SBarry Smith   IS               to,from;
21447c6ae99SBarry Smith   PetscErrorCode   ierr;
21547c6ae99SBarry Smith 
21647c6ae99SBarry Smith   PetscFunctionBegin;
217bff4a2f0SMatthew G. Knepley   if (stencil_type == DMDA_STENCIL_BOX && (bx == DM_BOUNDARY_MIRROR || by == DM_BOUNDARY_MIRROR)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Mirror boundary and box stencil");
21847c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
2193855c12bSBarry Smith #if !defined(PETSC_USE_64BIT_INDICES)
220bafee8b4SSatish Balay   if (((PetscInt64) M)*((PetscInt64) N)*((PetscInt64) dof) > (PetscInt64) PETSC_MPI_INT_MAX) SETERRQ3(comm,PETSC_ERR_INT_OVERFLOW,"Mesh of %D by %D by %D (dof) is too large for 32 bit indices",M,N,dof);
2213855c12bSBarry Smith #endif
22247c6ae99SBarry Smith 
22347c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
22447c6ae99SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
22547c6ae99SBarry Smith 
2267d310018SBarry Smith   dd->p = 1;
22747c6ae99SBarry Smith   if (m != PETSC_DECIDE) {
22847c6ae99SBarry Smith     if (m < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m);
22947c6ae99SBarry Smith     else if (m > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size);
23047c6ae99SBarry Smith   }
23147c6ae99SBarry Smith   if (n != PETSC_DECIDE) {
23247c6ae99SBarry Smith     if (n < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n);
23347c6ae99SBarry Smith     else if (n > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size);
23447c6ae99SBarry Smith   }
23547c6ae99SBarry Smith 
23647c6ae99SBarry Smith   if (m == PETSC_DECIDE || n == PETSC_DECIDE) {
23747c6ae99SBarry Smith     if (n != PETSC_DECIDE) {
23847c6ae99SBarry Smith       m = size/n;
23947c6ae99SBarry Smith     } else if (m != PETSC_DECIDE) {
24047c6ae99SBarry Smith       n = size/m;
24147c6ae99SBarry Smith     } else {
24247c6ae99SBarry Smith       /* try for squarish distribution */
243369cc0aeSBarry Smith       m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)N)));
24447c6ae99SBarry Smith       if (!m) m = 1;
24547c6ae99SBarry Smith       while (m > 0) {
24647c6ae99SBarry Smith         n = size/m;
24747c6ae99SBarry Smith         if (m*n == size) break;
24847c6ae99SBarry Smith         m--;
24947c6ae99SBarry Smith       }
25047c6ae99SBarry Smith       if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
25147c6ae99SBarry Smith     }
25247c6ae99SBarry Smith     if (m*n != size) SETERRQ(comm,PETSC_ERR_PLIB,"Unable to create partition, check the size of the communicator and input m and n ");
25347c6ae99SBarry Smith   } else if (m*n != size) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition");
25447c6ae99SBarry Smith 
25547c6ae99SBarry Smith   if (M < m) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m);
25647c6ae99SBarry Smith   if (N < n) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n);
25747c6ae99SBarry Smith 
25847c6ae99SBarry Smith   /*
25947c6ae99SBarry Smith      Determine locally owned region
26047c6ae99SBarry Smith      xs is the first local node number, x is the number of local nodes
26147c6ae99SBarry Smith   */
26247c6ae99SBarry Smith   if (!lx) {
263785e854fSJed Brown     ierr = PetscMalloc1(m, &dd->lx);CHKERRQ(ierr);
26447c6ae99SBarry Smith     lx   = dd->lx;
26547c6ae99SBarry Smith     for (i=0; i<m; i++) {
26647c6ae99SBarry Smith       lx[i] = M/m + ((M % m) > i);
26747c6ae99SBarry Smith     }
26847c6ae99SBarry Smith   }
26947c6ae99SBarry Smith   x  = lx[rank % m];
27047c6ae99SBarry Smith   xs = 0;
27147c6ae99SBarry Smith   for (i=0; i<(rank % m); i++) {
27247c6ae99SBarry Smith     xs += lx[i];
27347c6ae99SBarry Smith   }
274*76bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
27547c6ae99SBarry Smith     left = xs;
27647c6ae99SBarry Smith     for (i=(rank % m); i<m; i++) {
27747c6ae99SBarry Smith       left += lx[i];
27847c6ae99SBarry Smith     }
27947c6ae99SBarry Smith     if (left != M) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Sum of lx across processors not equal to M: %D %D",left,M);
280*76bd3646SJed Brown   }
28147c6ae99SBarry Smith 
28247c6ae99SBarry Smith   /*
28347c6ae99SBarry Smith      Determine locally owned region
28447c6ae99SBarry Smith      ys is the first local node number, y is the number of local nodes
28547c6ae99SBarry Smith   */
28647c6ae99SBarry Smith   if (!ly) {
287785e854fSJed Brown     ierr = PetscMalloc1(n, &dd->ly);CHKERRQ(ierr);
28847c6ae99SBarry Smith     ly   = dd->ly;
28947c6ae99SBarry Smith     for (i=0; i<n; i++) {
29047c6ae99SBarry Smith       ly[i] = N/n + ((N % n) > i);
29147c6ae99SBarry Smith     }
29247c6ae99SBarry Smith   }
29347c6ae99SBarry Smith   y  = ly[rank/m];
29447c6ae99SBarry Smith   ys = 0;
29547c6ae99SBarry Smith   for (i=0; i<(rank/m); i++) {
29647c6ae99SBarry Smith     ys += ly[i];
29747c6ae99SBarry Smith   }
298*76bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
29947c6ae99SBarry Smith     left = ys;
30047c6ae99SBarry Smith     for (i=(rank/m); i<n; i++) {
30147c6ae99SBarry Smith       left += ly[i];
30247c6ae99SBarry Smith     }
30347c6ae99SBarry Smith     if (left != N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Sum of ly across processors not equal to N: %D %D",left,N);
304*76bd3646SJed Brown   }
30547c6ae99SBarry Smith 
306bcea557cSEthan Coon   /*
307bcea557cSEthan Coon    check if the scatter requires more than one process neighbor or wraps around
308bcea557cSEthan Coon    the domain more than once
309bcea557cSEthan Coon   */
310bff4a2f0SMatthew G. Knepley   if ((x < s) && ((m > 1) || (bx == DM_BOUNDARY_PERIODIC))) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local x-width of domain x %D is smaller than stencil width s %D",x,s);
311bff4a2f0SMatthew G. Knepley   if ((y < s) && ((n > 1) || (by == DM_BOUNDARY_PERIODIC))) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local y-width of domain y %D is smaller than stencil width s %D",y,s);
31247c6ae99SBarry Smith   xe = xs + x;
31347c6ae99SBarry Smith   ye = ys + y;
31447c6ae99SBarry Smith 
315ce00eea3SSatish Balay   /* determine ghost region (Xs) and region scattered into (IXs)  */
316d9c9ebe5SPeter Brune   if (xs-s > 0) {
317d9c9ebe5SPeter Brune     Xs = xs - s; IXs = xs - s;
31888661749SPeter Brune   } else {
31988661749SPeter Brune     if (bx) {
32088661749SPeter Brune       Xs = xs - s;
32188661749SPeter Brune     } else {
32288661749SPeter Brune       Xs = 0;
32388661749SPeter Brune     }
32488661749SPeter Brune     IXs = 0;
32588661749SPeter Brune   }
326d9c9ebe5SPeter Brune   if (xe+s <= M) {
327d9c9ebe5SPeter Brune     Xe = xe + s; IXe = xe + s;
32888661749SPeter Brune   } else {
32988661749SPeter Brune     if (bx) {
330d9c9ebe5SPeter Brune       Xs = xs - s; Xe = xe + s;
33188661749SPeter Brune     } else {
33288661749SPeter Brune       Xe = M;
33388661749SPeter Brune     }
33488661749SPeter Brune     IXe = M;
33588661749SPeter Brune   }
33647c6ae99SBarry Smith 
337bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC || bx == DM_BOUNDARY_MIRROR) {
338d9c9ebe5SPeter Brune     IXs = xs - s;
339d9c9ebe5SPeter Brune     IXe = xe + s;
340d9c9ebe5SPeter Brune     Xs  = xs - s;
341d9c9ebe5SPeter Brune     Xe  = xe + s;
34288661749SPeter Brune   }
34347c6ae99SBarry Smith 
344d9c9ebe5SPeter Brune   if (ys-s > 0) {
345d9c9ebe5SPeter Brune     Ys = ys - s; IYs = ys - s;
34688661749SPeter Brune   } else {
34788661749SPeter Brune     if (by) {
34888661749SPeter Brune       Ys = ys - s;
34988661749SPeter Brune     } else {
35088661749SPeter Brune       Ys = 0;
35188661749SPeter Brune     }
35288661749SPeter Brune     IYs = 0;
35388661749SPeter Brune   }
354d9c9ebe5SPeter Brune   if (ye+s <= N) {
355d9c9ebe5SPeter Brune     Ye = ye + s; IYe = ye + s;
35688661749SPeter Brune   } else {
35788661749SPeter Brune     if (by) {
35888661749SPeter Brune       Ye = ye + s;
35988661749SPeter Brune     } else {
36088661749SPeter Brune       Ye = N;
36188661749SPeter Brune     }
36288661749SPeter Brune     IYe = N;
36388661749SPeter Brune   }
36488661749SPeter Brune 
365bff4a2f0SMatthew G. Knepley   if (by == DM_BOUNDARY_PERIODIC || by == DM_BOUNDARY_MIRROR) {
366d9c9ebe5SPeter Brune     IYs = ys - s;
367d9c9ebe5SPeter Brune     IYe = ye + s;
368d9c9ebe5SPeter Brune     Ys  = ys - s;
369d9c9ebe5SPeter Brune     Ye  = ye + s;
37088661749SPeter Brune   }
37188661749SPeter Brune 
37288661749SPeter Brune   /* stencil length in each direction */
373d9c9ebe5SPeter Brune   s_x = s;
374d9c9ebe5SPeter Brune   s_y = s;
37547c6ae99SBarry Smith 
37647c6ae99SBarry Smith   /* determine starting point of each processor */
37747c6ae99SBarry Smith   nn       = x*y;
378dcca6d9dSJed Brown   ierr     = PetscMalloc2(size+1,&bases,size,&ldims);CHKERRQ(ierr);
37947c6ae99SBarry Smith   ierr     = MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);CHKERRQ(ierr);
38047c6ae99SBarry Smith   bases[0] = 0;
38147c6ae99SBarry Smith   for (i=1; i<=size; i++) {
38247c6ae99SBarry Smith     bases[i] = ldims[i-1];
38347c6ae99SBarry Smith   }
38447c6ae99SBarry Smith   for (i=1; i<=size; i++) {
38547c6ae99SBarry Smith     bases[i] += bases[i-1];
38647c6ae99SBarry Smith   }
387ce00eea3SSatish Balay   base = bases[rank]*dof;
38847c6ae99SBarry Smith 
38947c6ae99SBarry Smith   /* allocate the base parallel and sequential vectors */
390ce00eea3SSatish Balay   dd->Nlocal = x*y*dof;
391b1fb7eb7SBarry Smith   ierr       = VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,NULL,&global);CHKERRQ(ierr);
392ce00eea3SSatish Balay   dd->nlocal = (Xe-Xs)*(Ye-Ys)*dof;
393b1fb7eb7SBarry Smith   ierr       = VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,NULL,&local);CHKERRQ(ierr);
39447c6ae99SBarry Smith 
3954104a7a0SPatrick Sanan   /* generate global to local vector scatter and local to global mapping*/
39647c6ae99SBarry Smith 
397ce00eea3SSatish Balay   /* global to local must include ghost points within the domain,
398ce00eea3SSatish Balay      but not ghost points outside the domain that aren't periodic */
3994104a7a0SPatrick Sanan   ierr  = PetscMalloc1((IXe-IXs)*(IYe-IYs),&idx);CHKERRQ(ierr);
400d9c9ebe5SPeter Brune   if (stencil_type == DMDA_STENCIL_BOX) {
401ce00eea3SSatish Balay     left  = IXs - Xs; right = left + (IXe-IXs);
402ce00eea3SSatish Balay     down  = IYs - Ys; up = down + (IYe-IYs);
403ce00eea3SSatish Balay     count = 0;
404ce00eea3SSatish Balay     for (i=down; i<up; i++) {
405ce00eea3SSatish Balay       for (j=left; j<right; j++) {
406ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
407ce00eea3SSatish Balay       }
408ce00eea3SSatish Balay     }
409ce00eea3SSatish Balay     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
410ce00eea3SSatish Balay 
41147c6ae99SBarry Smith   } else {
41247c6ae99SBarry Smith     /* must drop into cross shape region */
41347c6ae99SBarry Smith     /*       ---------|
41447c6ae99SBarry Smith             |  top    |
415ce00eea3SSatish Balay          |---         ---| up
41647c6ae99SBarry Smith          |   middle      |
41747c6ae99SBarry Smith          |               |
418ce00eea3SSatish Balay          ----         ---- down
41947c6ae99SBarry Smith             | bottom  |
42047c6ae99SBarry Smith             -----------
42147c6ae99SBarry Smith          Xs xs        xe Xe */
422ce00eea3SSatish Balay     left  = xs - Xs; right = left + x;
423ce00eea3SSatish Balay     down  = ys - Ys; up = down + y;
42447c6ae99SBarry Smith     count = 0;
425ce00eea3SSatish Balay     /* bottom */
426ce00eea3SSatish Balay     for (i=(IYs-Ys); i<down; i++) {
427ce00eea3SSatish Balay       for (j=left; j<right; j++) {
428ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
42947c6ae99SBarry Smith       }
43047c6ae99SBarry Smith     }
43147c6ae99SBarry Smith     /* middle */
43247c6ae99SBarry Smith     for (i=down; i<up; i++) {
433ce00eea3SSatish Balay       for (j=(IXs-Xs); j<(IXe-Xs); j++) {
434ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
43547c6ae99SBarry Smith       }
43647c6ae99SBarry Smith     }
43747c6ae99SBarry Smith     /* top */
438ce00eea3SSatish Balay     for (i=up; i<up+IYe-ye; i++) {
439ce00eea3SSatish Balay       for (j=left; j<right; j++) {
440ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
44147c6ae99SBarry Smith       }
44247c6ae99SBarry Smith     }
44347c6ae99SBarry Smith     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
44447c6ae99SBarry Smith   }
44547c6ae99SBarry Smith 
44647c6ae99SBarry Smith 
44747c6ae99SBarry Smith   /* determine who lies on each side of us stored in    n6 n7 n8
44847c6ae99SBarry Smith                                                         n3    n5
44947c6ae99SBarry Smith                                                         n0 n1 n2
45047c6ae99SBarry Smith   */
45147c6ae99SBarry Smith 
45247c6ae99SBarry Smith   /* Assume the Non-Periodic Case */
45347c6ae99SBarry Smith   n1 = rank - m;
45447c6ae99SBarry Smith   if (rank % m) {
45547c6ae99SBarry Smith     n0 = n1 - 1;
45647c6ae99SBarry Smith   } else {
45747c6ae99SBarry Smith     n0 = -1;
45847c6ae99SBarry Smith   }
45947c6ae99SBarry Smith   if ((rank+1) % m) {
46047c6ae99SBarry Smith     n2 = n1 + 1;
46147c6ae99SBarry Smith     n5 = rank + 1;
46247c6ae99SBarry Smith     n8 = rank + m + 1; if (n8 >= m*n) n8 = -1;
46347c6ae99SBarry Smith   } else {
46447c6ae99SBarry Smith     n2 = -1; n5 = -1; n8 = -1;
46547c6ae99SBarry Smith   }
46647c6ae99SBarry Smith   if (rank % m) {
46747c6ae99SBarry Smith     n3 = rank - 1;
46847c6ae99SBarry Smith     n6 = n3 + m; if (n6 >= m*n) n6 = -1;
46947c6ae99SBarry Smith   } else {
47047c6ae99SBarry Smith     n3 = -1; n6 = -1;
47147c6ae99SBarry Smith   }
47247c6ae99SBarry Smith   n7 = rank + m; if (n7 >= m*n) n7 = -1;
47347c6ae99SBarry Smith 
474bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC && by == DM_BOUNDARY_PERIODIC) {
47547c6ae99SBarry Smith     /* Modify for Periodic Cases */
47647c6ae99SBarry Smith     /* Handle all four corners */
47747c6ae99SBarry Smith     if ((n6 < 0) && (n7 < 0) && (n3 < 0)) n6 = m-1;
47847c6ae99SBarry Smith     if ((n8 < 0) && (n7 < 0) && (n5 < 0)) n8 = 0;
47947c6ae99SBarry Smith     if ((n2 < 0) && (n5 < 0) && (n1 < 0)) n2 = size-m;
48047c6ae99SBarry Smith     if ((n0 < 0) && (n3 < 0) && (n1 < 0)) n0 = size-1;
48147c6ae99SBarry Smith 
48247c6ae99SBarry Smith     /* Handle Top and Bottom Sides */
48347c6ae99SBarry Smith     if (n1 < 0) n1 = rank + m * (n-1);
48447c6ae99SBarry Smith     if (n7 < 0) n7 = rank - m * (n-1);
48547c6ae99SBarry Smith     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
48647c6ae99SBarry Smith     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
48747c6ae99SBarry Smith     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
48847c6ae99SBarry Smith     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
48947c6ae99SBarry Smith 
49047c6ae99SBarry Smith     /* Handle Left and Right Sides */
49147c6ae99SBarry Smith     if (n3 < 0) n3 = rank + (m-1);
49247c6ae99SBarry Smith     if (n5 < 0) n5 = rank - (m-1);
49347c6ae99SBarry Smith     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
49447c6ae99SBarry Smith     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
49547c6ae99SBarry Smith     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
49647c6ae99SBarry Smith     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
497bff4a2f0SMatthew G. Knepley   } else if (by == DM_BOUNDARY_PERIODIC) {  /* Handle Top and Bottom Sides */
498ce00eea3SSatish Balay     if (n1 < 0) n1 = rank + m * (n-1);
499ce00eea3SSatish Balay     if (n7 < 0) n7 = rank - m * (n-1);
500ce00eea3SSatish Balay     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
501ce00eea3SSatish Balay     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
502ce00eea3SSatish Balay     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
503ce00eea3SSatish Balay     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
504bff4a2f0SMatthew G. Knepley   } else if (bx == DM_BOUNDARY_PERIODIC) { /* Handle Left and Right Sides */
505ce00eea3SSatish Balay     if (n3 < 0) n3 = rank + (m-1);
506ce00eea3SSatish Balay     if (n5 < 0) n5 = rank - (m-1);
507ce00eea3SSatish Balay     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
508ce00eea3SSatish Balay     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
509ce00eea3SSatish Balay     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
510ce00eea3SSatish Balay     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
51147c6ae99SBarry Smith   }
512ce00eea3SSatish Balay 
513785e854fSJed Brown   ierr = PetscMalloc1(9,&dd->neighbors);CHKERRQ(ierr);
5148865f1eaSKarl Rupp 
51547c6ae99SBarry Smith   dd->neighbors[0] = n0;
51647c6ae99SBarry Smith   dd->neighbors[1] = n1;
51747c6ae99SBarry Smith   dd->neighbors[2] = n2;
51847c6ae99SBarry Smith   dd->neighbors[3] = n3;
51947c6ae99SBarry Smith   dd->neighbors[4] = rank;
52047c6ae99SBarry Smith   dd->neighbors[5] = n5;
52147c6ae99SBarry Smith   dd->neighbors[6] = n6;
52247c6ae99SBarry Smith   dd->neighbors[7] = n7;
52347c6ae99SBarry Smith   dd->neighbors[8] = n8;
52447c6ae99SBarry Smith 
525d9c9ebe5SPeter Brune   if (stencil_type == DMDA_STENCIL_STAR) {
52647c6ae99SBarry Smith     /* save corner processor numbers */
52747c6ae99SBarry Smith     sn0 = n0; sn2 = n2; sn6 = n6; sn8 = n8;
52847c6ae99SBarry Smith     n0  = n2 = n6 = n8 = -1;
52947c6ae99SBarry Smith   }
53047c6ae99SBarry Smith 
531785e854fSJed Brown   ierr = PetscMalloc1((Xe-Xs)*(Ye-Ys),&idx);CHKERRQ(ierr);
53247c6ae99SBarry Smith 
533ce00eea3SSatish Balay   nn = 0;
53447c6ae99SBarry Smith   xbase = bases[rank];
53547c6ae99SBarry Smith   for (i=1; i<=s_y; i++) {
53647c6ae99SBarry Smith     if (n0 >= 0) { /* left below */
537ce00eea3SSatish Balay       x_t = lx[n0 % m];
53847c6ae99SBarry Smith       y_t = ly[(n0/m)];
53947c6ae99SBarry Smith       s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
5408865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
54147c6ae99SBarry Smith     }
542ac119b13SBarry Smith 
54347c6ae99SBarry Smith     if (n1 >= 0) { /* directly below */
54447c6ae99SBarry Smith       x_t = x;
54547c6ae99SBarry Smith       y_t = ly[(n1/m)];
54647c6ae99SBarry Smith       s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
5478865f1eaSKarl Rupp       for (j=0; j<x_t; j++) idx[nn++] = s_t++;
548bff4a2f0SMatthew G. Knepley     } else if (by == DM_BOUNDARY_MIRROR) {
5498865f1eaSKarl Rupp       for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(s_y - i + 1)  + j;
55047c6ae99SBarry Smith     }
551ac119b13SBarry Smith 
55247c6ae99SBarry Smith     if (n2 >= 0) { /* right below */
553ce00eea3SSatish Balay       x_t = lx[n2 % m];
55447c6ae99SBarry Smith       y_t = ly[(n2/m)];
55547c6ae99SBarry Smith       s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
5568865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
55747c6ae99SBarry Smith     }
55847c6ae99SBarry Smith   }
55947c6ae99SBarry Smith 
56047c6ae99SBarry Smith   for (i=0; i<y; i++) {
56147c6ae99SBarry Smith     if (n3 >= 0) { /* directly left */
562ce00eea3SSatish Balay       x_t = lx[n3 % m];
56347c6ae99SBarry Smith       /* y_t = y; */
56447c6ae99SBarry Smith       s_t = bases[n3] + (i+1)*x_t - s_x;
5658865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
566bff4a2f0SMatthew G. Knepley     } else if (bx == DM_BOUNDARY_MIRROR) {
5678865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*i + s_x - j;
56847c6ae99SBarry Smith     }
56947c6ae99SBarry Smith 
5708865f1eaSKarl Rupp     for (j=0; j<x; j++) idx[nn++] = xbase++; /* interior */
57147c6ae99SBarry Smith 
57247c6ae99SBarry Smith     if (n5 >= 0) { /* directly right */
573ce00eea3SSatish Balay       x_t = lx[n5 % m];
57447c6ae99SBarry Smith       /* y_t = y; */
57547c6ae99SBarry Smith       s_t = bases[n5] + (i)*x_t;
5768865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
577bff4a2f0SMatthew G. Knepley     } else if (bx == DM_BOUNDARY_MIRROR) {
5788865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*(i + 1) - 2 - j;
57947c6ae99SBarry Smith     }
58047c6ae99SBarry Smith   }
58147c6ae99SBarry Smith 
58247c6ae99SBarry Smith   for (i=1; i<=s_y; i++) {
58347c6ae99SBarry Smith     if (n6 >= 0) { /* left above */
584ce00eea3SSatish Balay       x_t = lx[n6 % m];
58547c6ae99SBarry Smith       /* y_t = ly[(n6/m)]; */
58647c6ae99SBarry Smith       s_t = bases[n6] + (i)*x_t - s_x;
5878865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
58847c6ae99SBarry Smith     }
589ac119b13SBarry Smith 
59047c6ae99SBarry Smith     if (n7 >= 0) { /* directly above */
59147c6ae99SBarry Smith       x_t = x;
59247c6ae99SBarry Smith       /* y_t = ly[(n7/m)]; */
59347c6ae99SBarry Smith       s_t = bases[n7] + (i-1)*x_t;
5948865f1eaSKarl Rupp       for (j=0; j<x_t; j++) idx[nn++] = s_t++;
595bff4a2f0SMatthew G. Knepley     } else if (by == DM_BOUNDARY_MIRROR) {
5968865f1eaSKarl Rupp       for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(y - i - 1)  + j;
59747c6ae99SBarry Smith     }
598ac119b13SBarry Smith 
59947c6ae99SBarry Smith     if (n8 >= 0) { /* right above */
600ce00eea3SSatish Balay       x_t = lx[n8 % m];
60147c6ae99SBarry Smith       /* y_t = ly[(n8/m)]; */
60247c6ae99SBarry Smith       s_t = bases[n8] + (i-1)*x_t;
6038865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
60447c6ae99SBarry Smith     }
60547c6ae99SBarry Smith   }
60647c6ae99SBarry Smith 
607b1fb7eb7SBarry Smith   ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_USE_POINTER,&from);CHKERRQ(ierr);
6089448b7f1SJunchao Zhang   ierr = VecScatterCreate(global,from,local,to,&gtol);CHKERRQ(ierr);
6093bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)gtol);CHKERRQ(ierr);
610fcfd50ebSBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
611fcfd50ebSBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
61247c6ae99SBarry Smith 
613d9c9ebe5SPeter Brune   if (stencil_type == DMDA_STENCIL_STAR) {
614ce00eea3SSatish Balay     n0 = sn0; n2 = sn2; n6 = sn6; n8 = sn8;
615ce00eea3SSatish Balay   }
616ce00eea3SSatish Balay 
617288e7d53SBarry Smith   if (((stencil_type == DMDA_STENCIL_STAR)  || (bx && bx != DM_BOUNDARY_PERIODIC) || (by && by != DM_BOUNDARY_PERIODIC))) {
61847c6ae99SBarry Smith     /*
61947c6ae99SBarry Smith         Recompute the local to global mappings, this time keeping the
620ce00eea3SSatish Balay       information about the cross corner processor numbers and any ghosted
621ce00eea3SSatish Balay       but not periodic indices.
62247c6ae99SBarry Smith     */
62347c6ae99SBarry Smith     nn    = 0;
62447c6ae99SBarry Smith     xbase = bases[rank];
62547c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
62647c6ae99SBarry Smith       if (n0 >= 0) { /* left below */
627ce00eea3SSatish Balay         x_t = lx[n0 % m];
62847c6ae99SBarry Smith         y_t = ly[(n0/m)];
62947c6ae99SBarry Smith         s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
6308865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
631ce00eea3SSatish Balay       } else if (xs-Xs > 0 && ys-Ys > 0) {
6328865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = -1;
63347c6ae99SBarry Smith       }
63447c6ae99SBarry Smith       if (n1 >= 0) { /* directly below */
63547c6ae99SBarry Smith         x_t = x;
63647c6ae99SBarry Smith         y_t = ly[(n1/m)];
63747c6ae99SBarry Smith         s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
6388865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
639ce00eea3SSatish Balay       } else if (ys-Ys > 0) {
640bff4a2f0SMatthew G. Knepley         if (by == DM_BOUNDARY_MIRROR) {
6418865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(s_y - i + 1)  + j;
642624904c4SBarry Smith         } else {
6438865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = -1;
64447c6ae99SBarry Smith         }
645624904c4SBarry Smith       }
64647c6ae99SBarry Smith       if (n2 >= 0) { /* right below */
647ce00eea3SSatish Balay         x_t = lx[n2 % m];
64847c6ae99SBarry Smith         y_t = ly[(n2/m)];
64947c6ae99SBarry Smith         s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
6508865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
651ce00eea3SSatish Balay       } else if (Xe-xe> 0 && ys-Ys > 0) {
6528865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = -1;
65347c6ae99SBarry Smith       }
65447c6ae99SBarry Smith     }
65547c6ae99SBarry Smith 
65647c6ae99SBarry Smith     for (i=0; i<y; i++) {
65747c6ae99SBarry Smith       if (n3 >= 0) { /* directly left */
658ce00eea3SSatish Balay         x_t = lx[n3 % m];
65947c6ae99SBarry Smith         /* y_t = y; */
66047c6ae99SBarry Smith         s_t = bases[n3] + (i+1)*x_t - s_x;
6618865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
662ce00eea3SSatish Balay       } else if (xs-Xs > 0) {
663bff4a2f0SMatthew G. Knepley         if (bx == DM_BOUNDARY_MIRROR) {
6648865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*i + s_x - j;
665624904c4SBarry Smith         } else {
6668865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
66747c6ae99SBarry Smith         }
668624904c4SBarry Smith       }
66947c6ae99SBarry Smith 
6708865f1eaSKarl Rupp       for (j=0; j<x; j++) idx[nn++] = xbase++; /* interior */
67147c6ae99SBarry Smith 
67247c6ae99SBarry Smith       if (n5 >= 0) { /* directly right */
673ce00eea3SSatish Balay         x_t = lx[n5 % m];
67447c6ae99SBarry Smith         /* y_t = y; */
67547c6ae99SBarry Smith         s_t = bases[n5] + (i)*x_t;
6768865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
677ce00eea3SSatish Balay       } else if (Xe-xe > 0) {
678bff4a2f0SMatthew G. Knepley         if (bx == DM_BOUNDARY_MIRROR) {
6798865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*(i + 1) - 2 - j;
680624904c4SBarry Smith         } else {
6818865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
68247c6ae99SBarry Smith         }
68347c6ae99SBarry Smith       }
684624904c4SBarry Smith     }
68547c6ae99SBarry Smith 
68647c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
68747c6ae99SBarry Smith       if (n6 >= 0) { /* left above */
688ce00eea3SSatish Balay         x_t = lx[n6 % m];
68947c6ae99SBarry Smith         /* y_t = ly[(n6/m)]; */
69047c6ae99SBarry Smith         s_t = bases[n6] + (i)*x_t - s_x;
6918865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
692ce00eea3SSatish Balay       } else if (xs-Xs > 0 && Ye-ye > 0) {
6938865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = -1;
69447c6ae99SBarry Smith       }
69547c6ae99SBarry Smith       if (n7 >= 0) { /* directly above */
69647c6ae99SBarry Smith         x_t = x;
69747c6ae99SBarry Smith         /* y_t = ly[(n7/m)]; */
69847c6ae99SBarry Smith         s_t = bases[n7] + (i-1)*x_t;
6998865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
700ce00eea3SSatish Balay       } else if (Ye-ye > 0) {
701bff4a2f0SMatthew G. Knepley         if (by == DM_BOUNDARY_MIRROR) {
7028865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(y - i - 1)  + j;
703624904c4SBarry Smith         } else {
7048865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = -1;
70547c6ae99SBarry Smith         }
706624904c4SBarry Smith       }
70747c6ae99SBarry Smith       if (n8 >= 0) { /* right above */
708ce00eea3SSatish Balay         x_t = lx[n8 % m];
70947c6ae99SBarry Smith         /* y_t = ly[(n8/m)]; */
71047c6ae99SBarry Smith         s_t = bases[n8] + (i-1)*x_t;
7118865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
712ce00eea3SSatish Balay       } else if (Xe-xe > 0 && Ye-ye > 0) {
7138865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = -1;
71447c6ae99SBarry Smith       }
71547c6ae99SBarry Smith     }
71647c6ae99SBarry Smith   }
717ce00eea3SSatish Balay   /*
718ce00eea3SSatish Balay      Set the local to global ordering in the global vector, this allows use
719ce00eea3SSatish Balay      of VecSetValuesLocal().
720ce00eea3SSatish Balay   */
72145b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingCreate(comm,dof,nn,idx,PETSC_OWN_POINTER,&da->ltogmap);CHKERRQ(ierr);
7223bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)da->ltogmap);CHKERRQ(ierr);
72347c6ae99SBarry Smith 
724ce00eea3SSatish Balay   ierr  = PetscFree2(bases,ldims);CHKERRQ(ierr);
72547c6ae99SBarry Smith   dd->m = m;  dd->n  = n;
726ce00eea3SSatish Balay   /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */
727ce00eea3SSatish Balay   dd->xs = xs*dof; dd->xe = xe*dof; dd->ys = ys; dd->ye = ye; dd->zs = 0; dd->ze = 1;
728ce00eea3SSatish Balay   dd->Xs = Xs*dof; dd->Xe = Xe*dof; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = 0; dd->Ze = 1;
72947c6ae99SBarry Smith 
730fcfd50ebSBarry Smith   ierr = VecDestroy(&local);CHKERRQ(ierr);
731fcfd50ebSBarry Smith   ierr = VecDestroy(&global);CHKERRQ(ierr);
73247c6ae99SBarry Smith 
73347c6ae99SBarry Smith   dd->gtol      = gtol;
73447c6ae99SBarry Smith   dd->base      = base;
7359a42bb27SBarry Smith   da->ops->view = DMView_DA_2d;
7360298fd71SBarry Smith   dd->ltol      = NULL;
7370298fd71SBarry Smith   dd->ao        = NULL;
73847c6ae99SBarry Smith   PetscFunctionReturn(0);
73947c6ae99SBarry Smith }
74047c6ae99SBarry Smith 
74147c6ae99SBarry Smith /*@C
742aa219208SBarry Smith    DMDACreate2d -  Creates an object that will manage the communication of  two-dimensional
74347c6ae99SBarry Smith    regular array data that is distributed across some processors.
74447c6ae99SBarry Smith 
745d083f849SBarry Smith    Collective
74647c6ae99SBarry Smith 
74747c6ae99SBarry Smith    Input Parameters:
74847c6ae99SBarry Smith +  comm - MPI communicator
7491321219cSEthan Coon .  bx,by - type of ghost nodes the array have.
750bff4a2f0SMatthew G. Knepley          Use one of DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_PERIODIC.
751aa219208SBarry Smith .  stencil_type - stencil type.  Use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
752897f7067SBarry Smith .  M,N - global dimension in each direction of the array
75347c6ae99SBarry Smith .  m,n - corresponding number of processors in each dimension
75447c6ae99SBarry Smith          (or PETSC_DECIDE to have calculated)
75547c6ae99SBarry Smith .  dof - number of degrees of freedom per node
75647c6ae99SBarry Smith .  s - stencil width
75747c6ae99SBarry Smith -  lx, ly - arrays containing the number of nodes in each cell along
7580298fd71SBarry Smith            the x and y coordinates, or NULL. If non-null, these
75947c6ae99SBarry Smith            must be of length as m and n, and the corresponding
76047c6ae99SBarry Smith            m and n cannot be PETSC_DECIDE. The sum of the lx[] entries
76147c6ae99SBarry Smith            must be M, and the sum of the ly[] entries must be N.
76247c6ae99SBarry Smith 
76347c6ae99SBarry Smith    Output Parameter:
76447c6ae99SBarry Smith .  da - the resulting distributed array object
76547c6ae99SBarry Smith 
76647c6ae99SBarry Smith    Options Database Key:
767706a11cbSBarry Smith +  -dm_view - Calls DMView() at the conclusion of DMDACreate2d()
768e43dc8daSBarry Smith .  -da_grid_x <nx> - number of grid points in x direction
769e43dc8daSBarry Smith .  -da_grid_y <ny> - number of grid points in y direction
77047c6ae99SBarry Smith .  -da_processors_x <nx> - number of processors in x direction
77147c6ae99SBarry Smith .  -da_processors_y <ny> - number of processors in y direction
772e0f5d30fSBarry Smith .  -da_refine_x <rx> - refinement ratio in x direction
773e0f5d30fSBarry Smith .  -da_refine_y <ry> - refinement ratio in y direction
774e43dc8daSBarry Smith -  -da_refine <n> - refine the DMDA n times before creating
775e0f5d30fSBarry Smith 
77647c6ae99SBarry Smith 
77747c6ae99SBarry Smith    Level: beginner
77847c6ae99SBarry Smith 
77947c6ae99SBarry Smith    Notes:
780aa219208SBarry Smith    The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the
781aa219208SBarry Smith    standard 5-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes
78247c6ae99SBarry Smith    the standard 9-pt stencil.
78347c6ae99SBarry Smith 
784aa219208SBarry Smith    The array data itself is NOT stored in the DMDA, it is stored in Vec objects;
785564755cdSBarry Smith    The appropriate vector objects can be obtained with calls to DMCreateGlobalVector()
786564755cdSBarry Smith    and DMCreateLocalVector() and calls to VecDuplicate() if more are needed.
78747c6ae99SBarry Smith 
788897f7067SBarry Smith    You must call DMSetUp() after this call before using this DM.
789897f7067SBarry Smith 
790897f7067SBarry Smith    If you wish to use the options database to change values in the DMDA call DMSetFromOptions() after this call
791897f7067SBarry Smith    but before DMSetUp().
792897f7067SBarry Smith 
793aa219208SBarry Smith .seealso: DMDestroy(), DMView(), DMDACreate1d(), DMDACreate3d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(),
79499f0b487SRichard Tran Mills           DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMLocalToLocalBegin(), DMLocalToLocalEnd(), DMDASetRefinementFactor(),
7953bd220d7SPatrick Sanan           DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges(),
7963bd220d7SPatrick Sanan           DMStagCreate2d()
79747c6ae99SBarry Smith 
79847c6ae99SBarry Smith @*/
799fe16a2e9SBarry Smith 
800bff4a2f0SMatthew G. Knepley PetscErrorCode  DMDACreate2d(MPI_Comm comm,DMBoundaryType bx,DMBoundaryType by,DMDAStencilType stencil_type,
8019a42bb27SBarry Smith                              PetscInt M,PetscInt N,PetscInt m,PetscInt n,PetscInt dof,PetscInt s,const PetscInt lx[],const PetscInt ly[],DM *da)
80247c6ae99SBarry Smith {
80347c6ae99SBarry Smith   PetscErrorCode ierr;
80447c6ae99SBarry Smith 
80547c6ae99SBarry Smith   PetscFunctionBegin;
806aa219208SBarry Smith   ierr = DMDACreate(comm, da);CHKERRQ(ierr);
807c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(*da, 2);CHKERRQ(ierr);
808aa219208SBarry Smith   ierr = DMDASetSizes(*da, M, N, 1);CHKERRQ(ierr);
809aa219208SBarry Smith   ierr = DMDASetNumProcs(*da, m, n, PETSC_DECIDE);CHKERRQ(ierr);
810bff4a2f0SMatthew G. Knepley   ierr = DMDASetBoundaryType(*da, bx, by, DM_BOUNDARY_NONE);CHKERRQ(ierr);
811aa219208SBarry Smith   ierr = DMDASetDof(*da, dof);CHKERRQ(ierr);
812aa219208SBarry Smith   ierr = DMDASetStencilType(*da, stencil_type);CHKERRQ(ierr);
813aa219208SBarry Smith   ierr = DMDASetStencilWidth(*da, s);CHKERRQ(ierr);
8140298fd71SBarry Smith   ierr = DMDASetOwnershipRanges(*da, lx, ly, NULL);CHKERRQ(ierr);
81547c6ae99SBarry Smith   PetscFunctionReturn(0);
81647c6ae99SBarry Smith }
817