xref: /petsc/src/dm/impls/da/da2.c (revision e43dc8dad8d49be2730571154cf1a47f5fb69935)
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       M is number of grid points
14647c6ae99SBarry Smith       m is number of processors
14747c6ae99SBarry Smith 
14847c6ae99SBarry Smith */
1497087cfbeSBarry Smith PetscErrorCode  DMDASplitComm2d(MPI_Comm comm,PetscInt M,PetscInt N,PetscInt sw,MPI_Comm *outcomm)
15047c6ae99SBarry Smith {
15147c6ae99SBarry Smith   PetscErrorCode ierr;
15247c6ae99SBarry Smith   PetscInt       m,n = 0,x = 0,y = 0;
15347c6ae99SBarry Smith   PetscMPIInt    size,csize,rank;
15447c6ae99SBarry Smith 
15547c6ae99SBarry Smith   PetscFunctionBegin;
15647c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
15747c6ae99SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
15847c6ae99SBarry Smith 
15947c6ae99SBarry Smith   csize = 4*size;
16047c6ae99SBarry Smith   do {
16147c6ae99SBarry Smith     if (csize % 4) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cannot split communicator of size %d tried %d %D %D",size,csize,x,y);
16247c6ae99SBarry Smith     csize = csize/4;
16347c6ae99SBarry Smith 
164369cc0aeSBarry Smith     m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)csize)/((PetscReal)N)));
16547c6ae99SBarry Smith     if (!m) m = 1;
16647c6ae99SBarry Smith     while (m > 0) {
16747c6ae99SBarry Smith       n = csize/m;
16847c6ae99SBarry Smith       if (m*n == csize) break;
16947c6ae99SBarry Smith       m--;
17047c6ae99SBarry Smith     }
17147c6ae99SBarry Smith     if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
17247c6ae99SBarry Smith 
17347c6ae99SBarry Smith     x = M/m + ((M % m) > ((csize-1) % m));
17447c6ae99SBarry Smith     y = (N + (csize-1)/m)/n;
17547c6ae99SBarry Smith   } while ((x < 4 || y < 4) && csize > 1);
17647c6ae99SBarry Smith   if (size != csize) {
17747c6ae99SBarry Smith     MPI_Group   entire_group,sub_group;
17847c6ae99SBarry Smith     PetscMPIInt i,*groupies;
17947c6ae99SBarry Smith 
18047c6ae99SBarry Smith     ierr = MPI_Comm_group(comm,&entire_group);CHKERRQ(ierr);
181785e854fSJed Brown     ierr = PetscMalloc1(csize,&groupies);CHKERRQ(ierr);
18247c6ae99SBarry Smith     for (i=0; i<csize; i++) {
18347c6ae99SBarry Smith       groupies[i] = (rank/csize)*csize + i;
18447c6ae99SBarry Smith     }
18547c6ae99SBarry Smith     ierr = MPI_Group_incl(entire_group,csize,groupies,&sub_group);CHKERRQ(ierr);
18647c6ae99SBarry Smith     ierr = PetscFree(groupies);CHKERRQ(ierr);
18747c6ae99SBarry Smith     ierr = MPI_Comm_create(comm,sub_group,outcomm);CHKERRQ(ierr);
18847c6ae99SBarry Smith     ierr = MPI_Group_free(&entire_group);CHKERRQ(ierr);
18947c6ae99SBarry Smith     ierr = MPI_Group_free(&sub_group);CHKERRQ(ierr);
190aa219208SBarry Smith     ierr = PetscInfo1(0,"DMDASplitComm2d:Creating redundant coarse problems of size %d\n",csize);CHKERRQ(ierr);
19147c6ae99SBarry Smith   } else {
19247c6ae99SBarry Smith     *outcomm = comm;
19347c6ae99SBarry Smith   }
19447c6ae99SBarry Smith   PetscFunctionReturn(0);
19547c6ae99SBarry Smith }
19647c6ae99SBarry Smith 
19747c6ae99SBarry Smith #if defined(new)
19847c6ae99SBarry Smith /*
199aa219208SBarry Smith   DMDAGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix where local
200aa219208SBarry Smith     function lives on a DMDA
20147c6ae99SBarry Smith 
20247c6ae99SBarry Smith         y ~= (F(u + ha) - F(u))/h,
20347c6ae99SBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
20447c6ae99SBarry Smith         u = current iterate
20547c6ae99SBarry Smith         h = difference interval
20647c6ae99SBarry Smith */
207aa219208SBarry Smith PetscErrorCode DMDAGetDiagonal_MFFD(DM da,Vec U,Vec a)
20847c6ae99SBarry Smith {
20947c6ae99SBarry Smith   PetscScalar    h,*aa,*ww,v;
21047c6ae99SBarry Smith   PetscReal      epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
21147c6ae99SBarry Smith   PetscErrorCode ierr;
21247c6ae99SBarry Smith   PetscInt       gI,nI;
21347c6ae99SBarry Smith   MatStencil     stencil;
214aa219208SBarry Smith   DMDALocalInfo  info;
21547c6ae99SBarry Smith 
21647c6ae99SBarry Smith   PetscFunctionBegin;
21747c6ae99SBarry Smith   ierr = (*ctx->func)(0,U,a,ctx->funcctx);CHKERRQ(ierr);
21847c6ae99SBarry Smith   ierr = (*ctx->funcisetbase)(U,ctx->funcctx);CHKERRQ(ierr);
21947c6ae99SBarry Smith 
22047c6ae99SBarry Smith   ierr = VecGetArray(U,&ww);CHKERRQ(ierr);
22147c6ae99SBarry Smith   ierr = VecGetArray(a,&aa);CHKERRQ(ierr);
22247c6ae99SBarry Smith 
22347c6ae99SBarry Smith   nI = 0;
22447c6ae99SBarry Smith   h  = ww[gI];
22547c6ae99SBarry Smith   if (h == 0.0) h = 1.0;
22647c6ae99SBarry Smith   if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin;
22747c6ae99SBarry Smith   else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
22847c6ae99SBarry Smith   h *= epsilon;
22947c6ae99SBarry Smith 
23047c6ae99SBarry Smith   ww[gI] += h;
23147c6ae99SBarry Smith   ierr    = (*ctx->funci)(i,w,&v,ctx->funcctx);CHKERRQ(ierr);
23247c6ae99SBarry Smith   aa[nI]  = (v - aa[nI])/h;
23347c6ae99SBarry Smith   ww[gI] -= h;
23447c6ae99SBarry Smith   nI++;
2358865f1eaSKarl Rupp 
23647c6ae99SBarry Smith   ierr = VecRestoreArray(U,&ww);CHKERRQ(ierr);
23747c6ae99SBarry Smith   ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr);
23847c6ae99SBarry Smith   PetscFunctionReturn(0);
23947c6ae99SBarry Smith }
24047c6ae99SBarry Smith #endif
24147c6ae99SBarry Smith 
2427087cfbeSBarry Smith PetscErrorCode  DMSetUp_DA_2D(DM da)
24347c6ae99SBarry Smith {
24447c6ae99SBarry Smith   DM_DA            *dd = (DM_DA*)da->data;
24547c6ae99SBarry Smith   const PetscInt   M            = dd->M;
24647c6ae99SBarry Smith   const PetscInt   N            = dd->N;
24747c6ae99SBarry Smith   PetscInt         m            = dd->m;
24847c6ae99SBarry Smith   PetscInt         n            = dd->n;
24947c6ae99SBarry Smith   const PetscInt   dof          = dd->w;
25047c6ae99SBarry Smith   const PetscInt   s            = dd->s;
251bff4a2f0SMatthew G. Knepley   DMBoundaryType   bx           = dd->bx;
252bff4a2f0SMatthew G. Knepley   DMBoundaryType   by           = dd->by;
25319fd82e9SBarry Smith   DMDAStencilType  stencil_type = dd->stencil_type;
25447c6ae99SBarry Smith   PetscInt         *lx          = dd->lx;
25547c6ae99SBarry Smith   PetscInt         *ly          = dd->ly;
25647c6ae99SBarry Smith   MPI_Comm         comm;
25747c6ae99SBarry Smith   PetscMPIInt      rank,size;
258bd1fc5aeSBarry Smith   PetscInt         xs,xe,ys,ye,x,y,Xs,Xe,Ys,Ye,IXs,IXe,IYs,IYe;
2598ea3bf28SBarry Smith   PetscInt         up,down,left,right,i,n0,n1,n2,n3,n5,n6,n7,n8,*idx,nn;
260db87c5ecSEthan Coon   PetscInt         xbase,*bases,*ldims,j,x_t,y_t,s_t,base,count;
26147c6ae99SBarry Smith   PetscInt         s_x,s_y; /* s proportionalized to w */
26247c6ae99SBarry Smith   PetscInt         sn0 = 0,sn2 = 0,sn6 = 0,sn8 = 0;
26347c6ae99SBarry Smith   Vec              local,global;
264bd1fc5aeSBarry Smith   VecScatter       gtol;
26545b6f7e9SBarry Smith   IS               to,from;
26647c6ae99SBarry Smith   PetscErrorCode   ierr;
26747c6ae99SBarry Smith 
26847c6ae99SBarry Smith   PetscFunctionBegin;
269bff4a2f0SMatthew 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");
27047c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
2713855c12bSBarry Smith #if !defined(PETSC_USE_64BIT_INDICES)
272bafee8b4SSatish 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);
2733855c12bSBarry Smith #endif
27447c6ae99SBarry Smith 
27547c6ae99SBarry Smith   if (dof < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have 1 or more degrees of freedom per node: %D",dof);
27647c6ae99SBarry Smith   if (s < 0) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Stencil width cannot be negative: %D",s);
27747c6ae99SBarry Smith 
27847c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
27947c6ae99SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
28047c6ae99SBarry Smith 
2817d310018SBarry Smith   dd->p = 1;
28247c6ae99SBarry Smith   if (m != PETSC_DECIDE) {
28347c6ae99SBarry Smith     if (m < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m);
28447c6ae99SBarry Smith     else if (m > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size);
28547c6ae99SBarry Smith   }
28647c6ae99SBarry Smith   if (n != PETSC_DECIDE) {
28747c6ae99SBarry Smith     if (n < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n);
28847c6ae99SBarry Smith     else if (n > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size);
28947c6ae99SBarry Smith   }
29047c6ae99SBarry Smith 
29147c6ae99SBarry Smith   if (m == PETSC_DECIDE || n == PETSC_DECIDE) {
29247c6ae99SBarry Smith     if (n != PETSC_DECIDE) {
29347c6ae99SBarry Smith       m = size/n;
29447c6ae99SBarry Smith     } else if (m != PETSC_DECIDE) {
29547c6ae99SBarry Smith       n = size/m;
29647c6ae99SBarry Smith     } else {
29747c6ae99SBarry Smith       /* try for squarish distribution */
298369cc0aeSBarry Smith       m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)N)));
29947c6ae99SBarry Smith       if (!m) m = 1;
30047c6ae99SBarry Smith       while (m > 0) {
30147c6ae99SBarry Smith         n = size/m;
30247c6ae99SBarry Smith         if (m*n == size) break;
30347c6ae99SBarry Smith         m--;
30447c6ae99SBarry Smith       }
30547c6ae99SBarry Smith       if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
30647c6ae99SBarry Smith     }
30747c6ae99SBarry 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 ");
30847c6ae99SBarry Smith   } else if (m*n != size) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition");
30947c6ae99SBarry Smith 
31047c6ae99SBarry Smith   if (M < m) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m);
31147c6ae99SBarry Smith   if (N < n) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n);
31247c6ae99SBarry Smith 
31347c6ae99SBarry Smith   /*
31447c6ae99SBarry Smith      Determine locally owned region
31547c6ae99SBarry Smith      xs is the first local node number, x is the number of local nodes
31647c6ae99SBarry Smith   */
31747c6ae99SBarry Smith   if (!lx) {
318785e854fSJed Brown     ierr = PetscMalloc1(m, &dd->lx);CHKERRQ(ierr);
31947c6ae99SBarry Smith     lx   = dd->lx;
32047c6ae99SBarry Smith     for (i=0; i<m; i++) {
32147c6ae99SBarry Smith       lx[i] = M/m + ((M % m) > i);
32247c6ae99SBarry Smith     }
32347c6ae99SBarry Smith   }
32447c6ae99SBarry Smith   x  = lx[rank % m];
32547c6ae99SBarry Smith   xs = 0;
32647c6ae99SBarry Smith   for (i=0; i<(rank % m); i++) {
32747c6ae99SBarry Smith     xs += lx[i];
32847c6ae99SBarry Smith   }
32947c6ae99SBarry Smith #if defined(PETSC_USE_DEBUG)
33047c6ae99SBarry Smith   left = xs;
33147c6ae99SBarry Smith   for (i=(rank % m); i<m; i++) {
33247c6ae99SBarry Smith     left += lx[i];
33347c6ae99SBarry Smith   }
33447c6ae99SBarry 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);
33547c6ae99SBarry Smith #endif
33647c6ae99SBarry Smith 
33747c6ae99SBarry Smith   /*
33847c6ae99SBarry Smith      Determine locally owned region
33947c6ae99SBarry Smith      ys is the first local node number, y is the number of local nodes
34047c6ae99SBarry Smith   */
34147c6ae99SBarry Smith   if (!ly) {
342785e854fSJed Brown     ierr = PetscMalloc1(n, &dd->ly);CHKERRQ(ierr);
34347c6ae99SBarry Smith     ly   = dd->ly;
34447c6ae99SBarry Smith     for (i=0; i<n; i++) {
34547c6ae99SBarry Smith       ly[i] = N/n + ((N % n) > i);
34647c6ae99SBarry Smith     }
34747c6ae99SBarry Smith   }
34847c6ae99SBarry Smith   y  = ly[rank/m];
34947c6ae99SBarry Smith   ys = 0;
35047c6ae99SBarry Smith   for (i=0; i<(rank/m); i++) {
35147c6ae99SBarry Smith     ys += ly[i];
35247c6ae99SBarry Smith   }
35347c6ae99SBarry Smith #if defined(PETSC_USE_DEBUG)
35447c6ae99SBarry Smith   left = ys;
35547c6ae99SBarry Smith   for (i=(rank/m); i<n; i++) {
35647c6ae99SBarry Smith     left += ly[i];
35747c6ae99SBarry Smith   }
35847c6ae99SBarry 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);
35947c6ae99SBarry Smith #endif
36047c6ae99SBarry Smith 
361bcea557cSEthan Coon   /*
362bcea557cSEthan Coon    check if the scatter requires more than one process neighbor or wraps around
363bcea557cSEthan Coon    the domain more than once
364bcea557cSEthan Coon   */
365bff4a2f0SMatthew 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);
366bff4a2f0SMatthew 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);
36747c6ae99SBarry Smith   xe = xs + x;
36847c6ae99SBarry Smith   ye = ys + y;
36947c6ae99SBarry Smith 
370ce00eea3SSatish Balay   /* determine ghost region (Xs) and region scattered into (IXs)  */
371d9c9ebe5SPeter Brune   if (xs-s > 0) {
372d9c9ebe5SPeter Brune     Xs = xs - s; IXs = xs - s;
37388661749SPeter Brune   } else {
37488661749SPeter Brune     if (bx) {
37588661749SPeter Brune       Xs = xs - s;
37688661749SPeter Brune     } else {
37788661749SPeter Brune       Xs = 0;
37888661749SPeter Brune     }
37988661749SPeter Brune     IXs = 0;
38088661749SPeter Brune   }
381d9c9ebe5SPeter Brune   if (xe+s <= M) {
382d9c9ebe5SPeter Brune     Xe = xe + s; IXe = xe + s;
38388661749SPeter Brune   } else {
38488661749SPeter Brune     if (bx) {
385d9c9ebe5SPeter Brune       Xs = xs - s; Xe = xe + s;
38688661749SPeter Brune     } else {
38788661749SPeter Brune       Xe = M;
38888661749SPeter Brune     }
38988661749SPeter Brune     IXe = M;
39088661749SPeter Brune   }
39147c6ae99SBarry Smith 
392bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC || bx == DM_BOUNDARY_MIRROR) {
393d9c9ebe5SPeter Brune     IXs = xs - s;
394d9c9ebe5SPeter Brune     IXe = xe + s;
395d9c9ebe5SPeter Brune     Xs  = xs - s;
396d9c9ebe5SPeter Brune     Xe  = xe + s;
39788661749SPeter Brune   }
39847c6ae99SBarry Smith 
399d9c9ebe5SPeter Brune   if (ys-s > 0) {
400d9c9ebe5SPeter Brune     Ys = ys - s; IYs = ys - s;
40188661749SPeter Brune   } else {
40288661749SPeter Brune     if (by) {
40388661749SPeter Brune       Ys = ys - s;
40488661749SPeter Brune     } else {
40588661749SPeter Brune       Ys = 0;
40688661749SPeter Brune     }
40788661749SPeter Brune     IYs = 0;
40888661749SPeter Brune   }
409d9c9ebe5SPeter Brune   if (ye+s <= N) {
410d9c9ebe5SPeter Brune     Ye = ye + s; IYe = ye + s;
41188661749SPeter Brune   } else {
41288661749SPeter Brune     if (by) {
41388661749SPeter Brune       Ye = ye + s;
41488661749SPeter Brune     } else {
41588661749SPeter Brune       Ye = N;
41688661749SPeter Brune     }
41788661749SPeter Brune     IYe = N;
41888661749SPeter Brune   }
41988661749SPeter Brune 
420bff4a2f0SMatthew G. Knepley   if (by == DM_BOUNDARY_PERIODIC || by == DM_BOUNDARY_MIRROR) {
421d9c9ebe5SPeter Brune     IYs = ys - s;
422d9c9ebe5SPeter Brune     IYe = ye + s;
423d9c9ebe5SPeter Brune     Ys  = ys - s;
424d9c9ebe5SPeter Brune     Ye  = ye + s;
42588661749SPeter Brune   }
42688661749SPeter Brune 
42788661749SPeter Brune   /* stencil length in each direction */
428d9c9ebe5SPeter Brune   s_x = s;
429d9c9ebe5SPeter Brune   s_y = s;
43047c6ae99SBarry Smith 
43147c6ae99SBarry Smith   /* determine starting point of each processor */
43247c6ae99SBarry Smith   nn       = x*y;
433dcca6d9dSJed Brown   ierr     = PetscMalloc2(size+1,&bases,size,&ldims);CHKERRQ(ierr);
43447c6ae99SBarry Smith   ierr     = MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);CHKERRQ(ierr);
43547c6ae99SBarry Smith   bases[0] = 0;
43647c6ae99SBarry Smith   for (i=1; i<=size; i++) {
43747c6ae99SBarry Smith     bases[i] = ldims[i-1];
43847c6ae99SBarry Smith   }
43947c6ae99SBarry Smith   for (i=1; i<=size; i++) {
44047c6ae99SBarry Smith     bases[i] += bases[i-1];
44147c6ae99SBarry Smith   }
442ce00eea3SSatish Balay   base = bases[rank]*dof;
44347c6ae99SBarry Smith 
44447c6ae99SBarry Smith   /* allocate the base parallel and sequential vectors */
445ce00eea3SSatish Balay   dd->Nlocal = x*y*dof;
446b1fb7eb7SBarry Smith   ierr       = VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,NULL,&global);CHKERRQ(ierr);
447ce00eea3SSatish Balay   dd->nlocal = (Xe-Xs)*(Ye-Ys)*dof;
448b1fb7eb7SBarry Smith   ierr       = VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,NULL,&local);CHKERRQ(ierr);
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 */
4544104a7a0SPatrick Sanan   ierr  = PetscMalloc1((IXe-IXs)*(IYe-IYs),&idx);CHKERRQ(ierr);
455d9c9ebe5SPeter Brune   if (stencil_type == DMDA_STENCIL_BOX) {
456ce00eea3SSatish Balay     left  = IXs - Xs; right = left + (IXe-IXs);
457ce00eea3SSatish Balay     down  = IYs - Ys; up = down + (IYe-IYs);
458ce00eea3SSatish Balay     count = 0;
459ce00eea3SSatish Balay     for (i=down; i<up; i++) {
460ce00eea3SSatish Balay       for (j=left; j<right; j++) {
461ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
462ce00eea3SSatish Balay       }
463ce00eea3SSatish Balay     }
464ce00eea3SSatish Balay     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
465ce00eea3SSatish Balay 
46647c6ae99SBarry Smith   } else {
46747c6ae99SBarry Smith     /* must drop into cross shape region */
46847c6ae99SBarry Smith     /*       ---------|
46947c6ae99SBarry Smith             |  top    |
470ce00eea3SSatish Balay          |---         ---| up
47147c6ae99SBarry Smith          |   middle      |
47247c6ae99SBarry Smith          |               |
473ce00eea3SSatish Balay          ----         ---- down
47447c6ae99SBarry Smith             | bottom  |
47547c6ae99SBarry Smith             -----------
47647c6ae99SBarry Smith          Xs xs        xe Xe */
477ce00eea3SSatish Balay     left  = xs - Xs; right = left + x;
478ce00eea3SSatish Balay     down  = ys - Ys; up = down + y;
47947c6ae99SBarry Smith     count = 0;
480ce00eea3SSatish Balay     /* bottom */
481ce00eea3SSatish Balay     for (i=(IYs-Ys); i<down; i++) {
482ce00eea3SSatish Balay       for (j=left; j<right; j++) {
483ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
48447c6ae99SBarry Smith       }
48547c6ae99SBarry Smith     }
48647c6ae99SBarry Smith     /* middle */
48747c6ae99SBarry Smith     for (i=down; i<up; i++) {
488ce00eea3SSatish Balay       for (j=(IXs-Xs); j<(IXe-Xs); j++) {
489ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
49047c6ae99SBarry Smith       }
49147c6ae99SBarry Smith     }
49247c6ae99SBarry Smith     /* top */
493ce00eea3SSatish Balay     for (i=up; i<up+IYe-ye; i++) {
494ce00eea3SSatish Balay       for (j=left; j<right; j++) {
495ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
49647c6ae99SBarry Smith       }
49747c6ae99SBarry Smith     }
49847c6ae99SBarry Smith     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
49947c6ae99SBarry Smith   }
50047c6ae99SBarry Smith 
50147c6ae99SBarry Smith 
50247c6ae99SBarry Smith   /* determine who lies on each side of us stored in    n6 n7 n8
50347c6ae99SBarry Smith                                                         n3    n5
50447c6ae99SBarry Smith                                                         n0 n1 n2
50547c6ae99SBarry Smith   */
50647c6ae99SBarry Smith 
50747c6ae99SBarry Smith   /* Assume the Non-Periodic Case */
50847c6ae99SBarry Smith   n1 = rank - m;
50947c6ae99SBarry Smith   if (rank % m) {
51047c6ae99SBarry Smith     n0 = n1 - 1;
51147c6ae99SBarry Smith   } else {
51247c6ae99SBarry Smith     n0 = -1;
51347c6ae99SBarry Smith   }
51447c6ae99SBarry Smith   if ((rank+1) % m) {
51547c6ae99SBarry Smith     n2 = n1 + 1;
51647c6ae99SBarry Smith     n5 = rank + 1;
51747c6ae99SBarry Smith     n8 = rank + m + 1; if (n8 >= m*n) n8 = -1;
51847c6ae99SBarry Smith   } else {
51947c6ae99SBarry Smith     n2 = -1; n5 = -1; n8 = -1;
52047c6ae99SBarry Smith   }
52147c6ae99SBarry Smith   if (rank % m) {
52247c6ae99SBarry Smith     n3 = rank - 1;
52347c6ae99SBarry Smith     n6 = n3 + m; if (n6 >= m*n) n6 = -1;
52447c6ae99SBarry Smith   } else {
52547c6ae99SBarry Smith     n3 = -1; n6 = -1;
52647c6ae99SBarry Smith   }
52747c6ae99SBarry Smith   n7 = rank + m; if (n7 >= m*n) n7 = -1;
52847c6ae99SBarry Smith 
529bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC && by == DM_BOUNDARY_PERIODIC) {
53047c6ae99SBarry Smith     /* Modify for Periodic Cases */
53147c6ae99SBarry Smith     /* Handle all four corners */
53247c6ae99SBarry Smith     if ((n6 < 0) && (n7 < 0) && (n3 < 0)) n6 = m-1;
53347c6ae99SBarry Smith     if ((n8 < 0) && (n7 < 0) && (n5 < 0)) n8 = 0;
53447c6ae99SBarry Smith     if ((n2 < 0) && (n5 < 0) && (n1 < 0)) n2 = size-m;
53547c6ae99SBarry Smith     if ((n0 < 0) && (n3 < 0) && (n1 < 0)) n0 = size-1;
53647c6ae99SBarry Smith 
53747c6ae99SBarry Smith     /* Handle Top and Bottom Sides */
53847c6ae99SBarry Smith     if (n1 < 0) n1 = rank + m * (n-1);
53947c6ae99SBarry Smith     if (n7 < 0) n7 = rank - m * (n-1);
54047c6ae99SBarry Smith     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
54147c6ae99SBarry Smith     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
54247c6ae99SBarry Smith     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
54347c6ae99SBarry Smith     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
54447c6ae99SBarry Smith 
54547c6ae99SBarry Smith     /* Handle Left and Right Sides */
54647c6ae99SBarry Smith     if (n3 < 0) n3 = rank + (m-1);
54747c6ae99SBarry Smith     if (n5 < 0) n5 = rank - (m-1);
54847c6ae99SBarry Smith     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
54947c6ae99SBarry Smith     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
55047c6ae99SBarry Smith     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
55147c6ae99SBarry Smith     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
552bff4a2f0SMatthew G. Knepley   } else if (by == DM_BOUNDARY_PERIODIC) {  /* Handle Top and Bottom Sides */
553ce00eea3SSatish Balay     if (n1 < 0) n1 = rank + m * (n-1);
554ce00eea3SSatish Balay     if (n7 < 0) n7 = rank - m * (n-1);
555ce00eea3SSatish Balay     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
556ce00eea3SSatish Balay     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
557ce00eea3SSatish Balay     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
558ce00eea3SSatish Balay     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
559bff4a2f0SMatthew G. Knepley   } else if (bx == DM_BOUNDARY_PERIODIC) { /* Handle Left and Right Sides */
560ce00eea3SSatish Balay     if (n3 < 0) n3 = rank + (m-1);
561ce00eea3SSatish Balay     if (n5 < 0) n5 = rank - (m-1);
562ce00eea3SSatish Balay     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
563ce00eea3SSatish Balay     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
564ce00eea3SSatish Balay     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
565ce00eea3SSatish Balay     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
56647c6ae99SBarry Smith   }
567ce00eea3SSatish Balay 
568785e854fSJed Brown   ierr = PetscMalloc1(9,&dd->neighbors);CHKERRQ(ierr);
5698865f1eaSKarl Rupp 
57047c6ae99SBarry Smith   dd->neighbors[0] = n0;
57147c6ae99SBarry Smith   dd->neighbors[1] = n1;
57247c6ae99SBarry Smith   dd->neighbors[2] = n2;
57347c6ae99SBarry Smith   dd->neighbors[3] = n3;
57447c6ae99SBarry Smith   dd->neighbors[4] = rank;
57547c6ae99SBarry Smith   dd->neighbors[5] = n5;
57647c6ae99SBarry Smith   dd->neighbors[6] = n6;
57747c6ae99SBarry Smith   dd->neighbors[7] = n7;
57847c6ae99SBarry Smith   dd->neighbors[8] = n8;
57947c6ae99SBarry Smith 
580d9c9ebe5SPeter Brune   if (stencil_type == DMDA_STENCIL_STAR) {
58147c6ae99SBarry Smith     /* save corner processor numbers */
58247c6ae99SBarry Smith     sn0 = n0; sn2 = n2; sn6 = n6; sn8 = n8;
58347c6ae99SBarry Smith     n0  = n2 = n6 = n8 = -1;
58447c6ae99SBarry Smith   }
58547c6ae99SBarry Smith 
586785e854fSJed Brown   ierr = PetscMalloc1((Xe-Xs)*(Ye-Ys),&idx);CHKERRQ(ierr);
58747c6ae99SBarry Smith 
588ce00eea3SSatish Balay   nn = 0;
58947c6ae99SBarry Smith   xbase = bases[rank];
59047c6ae99SBarry Smith   for (i=1; i<=s_y; i++) {
59147c6ae99SBarry Smith     if (n0 >= 0) { /* left below */
592ce00eea3SSatish Balay       x_t = lx[n0 % m];
59347c6ae99SBarry Smith       y_t = ly[(n0/m)];
59447c6ae99SBarry Smith       s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
5958865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
59647c6ae99SBarry Smith     }
597ac119b13SBarry Smith 
59847c6ae99SBarry Smith     if (n1 >= 0) { /* directly below */
59947c6ae99SBarry Smith       x_t = x;
60047c6ae99SBarry Smith       y_t = ly[(n1/m)];
60147c6ae99SBarry Smith       s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
6028865f1eaSKarl Rupp       for (j=0; j<x_t; j++) idx[nn++] = s_t++;
603bff4a2f0SMatthew G. Knepley     } else if (by == DM_BOUNDARY_MIRROR) {
6048865f1eaSKarl Rupp       for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(s_y - i + 1)  + j;
60547c6ae99SBarry Smith     }
606ac119b13SBarry Smith 
60747c6ae99SBarry Smith     if (n2 >= 0) { /* right below */
608ce00eea3SSatish Balay       x_t = lx[n2 % m];
60947c6ae99SBarry Smith       y_t = ly[(n2/m)];
61047c6ae99SBarry Smith       s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
6118865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
61247c6ae99SBarry Smith     }
61347c6ae99SBarry Smith   }
61447c6ae99SBarry Smith 
61547c6ae99SBarry Smith   for (i=0; i<y; i++) {
61647c6ae99SBarry Smith     if (n3 >= 0) { /* directly left */
617ce00eea3SSatish Balay       x_t = lx[n3 % m];
61847c6ae99SBarry Smith       /* y_t = y; */
61947c6ae99SBarry Smith       s_t = bases[n3] + (i+1)*x_t - s_x;
6208865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
621bff4a2f0SMatthew G. Knepley     } else if (bx == DM_BOUNDARY_MIRROR) {
6228865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*i + s_x - j;
62347c6ae99SBarry Smith     }
62447c6ae99SBarry Smith 
6258865f1eaSKarl Rupp     for (j=0; j<x; j++) idx[nn++] = xbase++; /* interior */
62647c6ae99SBarry Smith 
62747c6ae99SBarry Smith     if (n5 >= 0) { /* directly right */
628ce00eea3SSatish Balay       x_t = lx[n5 % m];
62947c6ae99SBarry Smith       /* y_t = y; */
63047c6ae99SBarry Smith       s_t = bases[n5] + (i)*x_t;
6318865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
632bff4a2f0SMatthew G. Knepley     } else if (bx == DM_BOUNDARY_MIRROR) {
6338865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*(i + 1) - 2 - j;
63447c6ae99SBarry Smith     }
63547c6ae99SBarry Smith   }
63647c6ae99SBarry Smith 
63747c6ae99SBarry Smith   for (i=1; i<=s_y; i++) {
63847c6ae99SBarry Smith     if (n6 >= 0) { /* left above */
639ce00eea3SSatish Balay       x_t = lx[n6 % m];
64047c6ae99SBarry Smith       /* y_t = ly[(n6/m)]; */
64147c6ae99SBarry Smith       s_t = bases[n6] + (i)*x_t - s_x;
6428865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
64347c6ae99SBarry Smith     }
644ac119b13SBarry Smith 
64547c6ae99SBarry Smith     if (n7 >= 0) { /* directly above */
64647c6ae99SBarry Smith       x_t = x;
64747c6ae99SBarry Smith       /* y_t = ly[(n7/m)]; */
64847c6ae99SBarry Smith       s_t = bases[n7] + (i-1)*x_t;
6498865f1eaSKarl Rupp       for (j=0; j<x_t; j++) idx[nn++] = s_t++;
650bff4a2f0SMatthew G. Knepley     } else if (by == DM_BOUNDARY_MIRROR) {
6518865f1eaSKarl Rupp       for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(y - i - 1)  + j;
65247c6ae99SBarry Smith     }
653ac119b13SBarry Smith 
65447c6ae99SBarry Smith     if (n8 >= 0) { /* right above */
655ce00eea3SSatish Balay       x_t = lx[n8 % m];
65647c6ae99SBarry Smith       /* y_t = ly[(n8/m)]; */
65747c6ae99SBarry Smith       s_t = bases[n8] + (i-1)*x_t;
6588865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
65947c6ae99SBarry Smith     }
66047c6ae99SBarry Smith   }
66147c6ae99SBarry Smith 
662b1fb7eb7SBarry Smith   ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_USE_POINTER,&from);CHKERRQ(ierr);
66347c6ae99SBarry Smith   ierr = VecScatterCreate(global,from,local,to,&gtol);CHKERRQ(ierr);
6643bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)gtol);CHKERRQ(ierr);
665fcfd50ebSBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
666fcfd50ebSBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
66747c6ae99SBarry Smith 
668d9c9ebe5SPeter Brune   if (stencil_type == DMDA_STENCIL_STAR) {
669ce00eea3SSatish Balay     n0 = sn0; n2 = sn2; n6 = sn6; n8 = sn8;
670ce00eea3SSatish Balay   }
671ce00eea3SSatish Balay 
672288e7d53SBarry Smith   if (((stencil_type == DMDA_STENCIL_STAR)  || (bx && bx != DM_BOUNDARY_PERIODIC) || (by && by != DM_BOUNDARY_PERIODIC))) {
67347c6ae99SBarry Smith     /*
67447c6ae99SBarry Smith         Recompute the local to global mappings, this time keeping the
675ce00eea3SSatish Balay       information about the cross corner processor numbers and any ghosted
676ce00eea3SSatish Balay       but not periodic indices.
67747c6ae99SBarry Smith     */
67847c6ae99SBarry Smith     nn    = 0;
67947c6ae99SBarry Smith     xbase = bases[rank];
68047c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
68147c6ae99SBarry Smith       if (n0 >= 0) { /* left below */
682ce00eea3SSatish Balay         x_t = lx[n0 % m];
68347c6ae99SBarry Smith         y_t = ly[(n0/m)];
68447c6ae99SBarry Smith         s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
6858865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
686ce00eea3SSatish Balay       } else if (xs-Xs > 0 && ys-Ys > 0) {
6878865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = -1;
68847c6ae99SBarry Smith       }
68947c6ae99SBarry Smith       if (n1 >= 0) { /* directly below */
69047c6ae99SBarry Smith         x_t = x;
69147c6ae99SBarry Smith         y_t = ly[(n1/m)];
69247c6ae99SBarry Smith         s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
6938865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
694ce00eea3SSatish Balay       } else if (ys-Ys > 0) {
695bff4a2f0SMatthew G. Knepley         if (by == DM_BOUNDARY_MIRROR) {
6968865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(s_y - i + 1)  + j;
697624904c4SBarry Smith         } else {
6988865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = -1;
69947c6ae99SBarry Smith         }
700624904c4SBarry Smith       }
70147c6ae99SBarry Smith       if (n2 >= 0) { /* right below */
702ce00eea3SSatish Balay         x_t = lx[n2 % m];
70347c6ae99SBarry Smith         y_t = ly[(n2/m)];
70447c6ae99SBarry Smith         s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
7058865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
706ce00eea3SSatish Balay       } else if (Xe-xe> 0 && ys-Ys > 0) {
7078865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = -1;
70847c6ae99SBarry Smith       }
70947c6ae99SBarry Smith     }
71047c6ae99SBarry Smith 
71147c6ae99SBarry Smith     for (i=0; i<y; i++) {
71247c6ae99SBarry Smith       if (n3 >= 0) { /* directly left */
713ce00eea3SSatish Balay         x_t = lx[n3 % m];
71447c6ae99SBarry Smith         /* y_t = y; */
71547c6ae99SBarry Smith         s_t = bases[n3] + (i+1)*x_t - s_x;
7168865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
717ce00eea3SSatish Balay       } else if (xs-Xs > 0) {
718bff4a2f0SMatthew G. Knepley         if (bx == DM_BOUNDARY_MIRROR) {
7198865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*i + s_x - j;
720624904c4SBarry Smith         } else {
7218865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
72247c6ae99SBarry Smith         }
723624904c4SBarry Smith       }
72447c6ae99SBarry Smith 
7258865f1eaSKarl Rupp       for (j=0; j<x; j++) idx[nn++] = xbase++; /* interior */
72647c6ae99SBarry Smith 
72747c6ae99SBarry Smith       if (n5 >= 0) { /* directly right */
728ce00eea3SSatish Balay         x_t = lx[n5 % m];
72947c6ae99SBarry Smith         /* y_t = y; */
73047c6ae99SBarry Smith         s_t = bases[n5] + (i)*x_t;
7318865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
732ce00eea3SSatish Balay       } else if (Xe-xe > 0) {
733bff4a2f0SMatthew G. Knepley         if (bx == DM_BOUNDARY_MIRROR) {
7348865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*(i + 1) - 2 - j;
735624904c4SBarry Smith         } else {
7368865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
73747c6ae99SBarry Smith         }
73847c6ae99SBarry Smith       }
739624904c4SBarry Smith     }
74047c6ae99SBarry Smith 
74147c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
74247c6ae99SBarry Smith       if (n6 >= 0) { /* left above */
743ce00eea3SSatish Balay         x_t = lx[n6 % m];
74447c6ae99SBarry Smith         /* y_t = ly[(n6/m)]; */
74547c6ae99SBarry Smith         s_t = bases[n6] + (i)*x_t - s_x;
7468865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
747ce00eea3SSatish Balay       } else if (xs-Xs > 0 && Ye-ye > 0) {
7488865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = -1;
74947c6ae99SBarry Smith       }
75047c6ae99SBarry Smith       if (n7 >= 0) { /* directly above */
75147c6ae99SBarry Smith         x_t = x;
75247c6ae99SBarry Smith         /* y_t = ly[(n7/m)]; */
75347c6ae99SBarry Smith         s_t = bases[n7] + (i-1)*x_t;
7548865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
755ce00eea3SSatish Balay       } else if (Ye-ye > 0) {
756bff4a2f0SMatthew G. Knepley         if (by == DM_BOUNDARY_MIRROR) {
7578865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(y - i - 1)  + j;
758624904c4SBarry Smith         } else {
7598865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = -1;
76047c6ae99SBarry Smith         }
761624904c4SBarry Smith       }
76247c6ae99SBarry Smith       if (n8 >= 0) { /* right above */
763ce00eea3SSatish Balay         x_t = lx[n8 % m];
76447c6ae99SBarry Smith         /* y_t = ly[(n8/m)]; */
76547c6ae99SBarry Smith         s_t = bases[n8] + (i-1)*x_t;
7668865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
767ce00eea3SSatish Balay       } else if (Xe-xe > 0 && Ye-ye > 0) {
7688865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = -1;
76947c6ae99SBarry Smith       }
77047c6ae99SBarry Smith     }
77147c6ae99SBarry Smith   }
772ce00eea3SSatish Balay   /*
773ce00eea3SSatish Balay      Set the local to global ordering in the global vector, this allows use
774ce00eea3SSatish Balay      of VecSetValuesLocal().
775ce00eea3SSatish Balay   */
77645b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingCreate(comm,dof,nn,idx,PETSC_OWN_POINTER,&da->ltogmap);CHKERRQ(ierr);
7773bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)da->ltogmap);CHKERRQ(ierr);
77847c6ae99SBarry Smith 
779ce00eea3SSatish Balay   ierr  = PetscFree2(bases,ldims);CHKERRQ(ierr);
78047c6ae99SBarry Smith   dd->m = m;  dd->n  = n;
781ce00eea3SSatish Balay   /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */
782ce00eea3SSatish Balay   dd->xs = xs*dof; dd->xe = xe*dof; dd->ys = ys; dd->ye = ye; dd->zs = 0; dd->ze = 1;
783ce00eea3SSatish Balay   dd->Xs = Xs*dof; dd->Xe = Xe*dof; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = 0; dd->Ze = 1;
78447c6ae99SBarry Smith 
785fcfd50ebSBarry Smith   ierr = VecDestroy(&local);CHKERRQ(ierr);
786fcfd50ebSBarry Smith   ierr = VecDestroy(&global);CHKERRQ(ierr);
78747c6ae99SBarry Smith 
78847c6ae99SBarry Smith   dd->gtol      = gtol;
78947c6ae99SBarry Smith   dd->base      = base;
7909a42bb27SBarry Smith   da->ops->view = DMView_DA_2d;
7910298fd71SBarry Smith   dd->ltol      = NULL;
7920298fd71SBarry Smith   dd->ao        = NULL;
79347c6ae99SBarry Smith   PetscFunctionReturn(0);
79447c6ae99SBarry Smith }
79547c6ae99SBarry Smith 
79647c6ae99SBarry Smith /*@C
797aa219208SBarry Smith    DMDACreate2d -  Creates an object that will manage the communication of  two-dimensional
79847c6ae99SBarry Smith    regular array data that is distributed across some processors.
79947c6ae99SBarry Smith 
80047c6ae99SBarry Smith    Collective on MPI_Comm
80147c6ae99SBarry Smith 
80247c6ae99SBarry Smith    Input Parameters:
80347c6ae99SBarry Smith +  comm - MPI communicator
8041321219cSEthan Coon .  bx,by - type of ghost nodes the array have.
805bff4a2f0SMatthew G. Knepley          Use one of DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_PERIODIC.
806aa219208SBarry Smith .  stencil_type - stencil type.  Use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
807897f7067SBarry Smith .  M,N - global dimension in each direction of the array
80847c6ae99SBarry Smith .  m,n - corresponding number of processors in each dimension
80947c6ae99SBarry Smith          (or PETSC_DECIDE to have calculated)
81047c6ae99SBarry Smith .  dof - number of degrees of freedom per node
81147c6ae99SBarry Smith .  s - stencil width
81247c6ae99SBarry Smith -  lx, ly - arrays containing the number of nodes in each cell along
8130298fd71SBarry Smith            the x and y coordinates, or NULL. If non-null, these
81447c6ae99SBarry Smith            must be of length as m and n, and the corresponding
81547c6ae99SBarry Smith            m and n cannot be PETSC_DECIDE. The sum of the lx[] entries
81647c6ae99SBarry Smith            must be M, and the sum of the ly[] entries must be N.
81747c6ae99SBarry Smith 
81847c6ae99SBarry Smith    Output Parameter:
81947c6ae99SBarry Smith .  da - the resulting distributed array object
82047c6ae99SBarry Smith 
82147c6ae99SBarry Smith    Options Database Key:
822706a11cbSBarry Smith +  -dm_view - Calls DMView() at the conclusion of DMDACreate2d()
823*e43dc8daSBarry Smith .  -da_grid_x <nx> - number of grid points in x direction
824*e43dc8daSBarry Smith .  -da_grid_y <ny> - number of grid points in y direction
82547c6ae99SBarry Smith .  -da_processors_x <nx> - number of processors in x direction
82647c6ae99SBarry Smith .  -da_processors_y <ny> - number of processors in y direction
827e0f5d30fSBarry Smith .  -da_refine_x <rx> - refinement ratio in x direction
828e0f5d30fSBarry Smith .  -da_refine_y <ry> - refinement ratio in y direction
829*e43dc8daSBarry Smith -  -da_refine <n> - refine the DMDA n times before creating
830e0f5d30fSBarry Smith 
83147c6ae99SBarry Smith 
83247c6ae99SBarry Smith    Level: beginner
83347c6ae99SBarry Smith 
83447c6ae99SBarry Smith    Notes:
835aa219208SBarry Smith    The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the
836aa219208SBarry Smith    standard 5-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes
83747c6ae99SBarry Smith    the standard 9-pt stencil.
83847c6ae99SBarry Smith 
839aa219208SBarry Smith    The array data itself is NOT stored in the DMDA, it is stored in Vec objects;
840564755cdSBarry Smith    The appropriate vector objects can be obtained with calls to DMCreateGlobalVector()
841564755cdSBarry Smith    and DMCreateLocalVector() and calls to VecDuplicate() if more are needed.
84247c6ae99SBarry Smith 
843897f7067SBarry Smith    You must call DMSetUp() after this call before using this DM.
844897f7067SBarry Smith 
845897f7067SBarry Smith    If you wish to use the options database to change values in the DMDA call DMSetFromOptions() after this call
846897f7067SBarry Smith    but before DMSetUp().
847897f7067SBarry Smith 
84847c6ae99SBarry Smith .keywords: distributed array, create, two-dimensional
84947c6ae99SBarry Smith 
850aa219208SBarry Smith .seealso: DMDestroy(), DMView(), DMDACreate1d(), DMDACreate3d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(),
85199f0b487SRichard Tran Mills           DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMLocalToLocalBegin(), DMLocalToLocalEnd(), DMDASetRefinementFactor(),
852d461ba97SBarry Smith           DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges()
85347c6ae99SBarry Smith 
85447c6ae99SBarry Smith @*/
855fe16a2e9SBarry Smith 
856bff4a2f0SMatthew G. Knepley PetscErrorCode  DMDACreate2d(MPI_Comm comm,DMBoundaryType bx,DMBoundaryType by,DMDAStencilType stencil_type,
8579a42bb27SBarry Smith                           PetscInt M,PetscInt N,PetscInt m,PetscInt n,PetscInt dof,PetscInt s,const PetscInt lx[],const PetscInt ly[],DM *da)
85847c6ae99SBarry Smith {
85947c6ae99SBarry Smith   PetscErrorCode ierr;
86047c6ae99SBarry Smith 
86147c6ae99SBarry Smith   PetscFunctionBegin;
862aa219208SBarry Smith   ierr = DMDACreate(comm, da);CHKERRQ(ierr);
863c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(*da, 2);CHKERRQ(ierr);
864aa219208SBarry Smith   ierr = DMDASetSizes(*da, M, N, 1);CHKERRQ(ierr);
865aa219208SBarry Smith   ierr = DMDASetNumProcs(*da, m, n, PETSC_DECIDE);CHKERRQ(ierr);
866bff4a2f0SMatthew G. Knepley   ierr = DMDASetBoundaryType(*da, bx, by, DM_BOUNDARY_NONE);CHKERRQ(ierr);
867aa219208SBarry Smith   ierr = DMDASetDof(*da, dof);CHKERRQ(ierr);
868aa219208SBarry Smith   ierr = DMDASetStencilType(*da, stencil_type);CHKERRQ(ierr);
869aa219208SBarry Smith   ierr = DMDASetStencilWidth(*da, s);CHKERRQ(ierr);
8700298fd71SBarry Smith   ierr = DMDASetOwnershipRanges(*da, lx, ly, NULL);CHKERRQ(ierr);
87147c6ae99SBarry Smith   PetscFunctionReturn(0);
87247c6ae99SBarry Smith }
873