19a42bb27SBarry Smith 24035e84dSBarry Smith #include <petsc-private/dmdaimpl.h> /*I "petscdmda.h" I*/ 39804daf3SBarry Smith #include <petscdraw.h> 447c6ae99SBarry Smith 547c6ae99SBarry Smith #undef __FUNCT__ 69a42bb27SBarry Smith #define __FUNCT__ "DMView_DA_2d" 79a42bb27SBarry Smith PetscErrorCode DMView_DA_2d(DM da,PetscViewer viewer) 847c6ae99SBarry Smith { 947c6ae99SBarry Smith PetscErrorCode ierr; 1047c6ae99SBarry Smith PetscMPIInt rank; 119a42bb27SBarry Smith PetscBool iascii,isdraw,isbinary; 1247c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 139a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 149a42bb27SBarry Smith PetscBool ismatlab; 159a42bb27SBarry Smith #endif 1647c6ae99SBarry Smith 1747c6ae99SBarry Smith PetscFunctionBegin; 18ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)da),&rank);CHKERRQ(ierr); 1947c6ae99SBarry Smith 20251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 21251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 22251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 239a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 24251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERMATLAB,&ismatlab);CHKERRQ(ierr); 259a42bb27SBarry Smith #endif 2647c6ae99SBarry Smith if (iascii) { 2747c6ae99SBarry Smith PetscViewerFormat format; 2847c6ae99SBarry Smith 2947c6ae99SBarry Smith ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr); 3047c6ae99SBarry Smith if (format != PETSC_VIEWER_ASCII_VTK && format != PETSC_VIEWER_ASCII_VTK_CELL) { 31aa219208SBarry Smith DMDALocalInfo info; 32aa219208SBarry Smith ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr); 337b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 3447c6ae99SBarry 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); 3547c6ae99SBarry 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); 3647c6ae99SBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 377b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr); 383da9ae13SJed Brown } else { 393da9ae13SJed Brown ierr = DMView_DA_VTK(da,viewer);CHKERRQ(ierr); 4047c6ae99SBarry Smith } 4147c6ae99SBarry Smith } else if (isdraw) { 4247c6ae99SBarry Smith PetscDraw draw; 4347c6ae99SBarry Smith double ymin = -1*dd->s-1,ymax = dd->N+dd->s; 4447c6ae99SBarry Smith double xmin = -1*dd->s-1,xmax = dd->M+dd->s; 4547c6ae99SBarry Smith double x,y; 468ea3bf28SBarry Smith PetscInt base; 478ea3bf28SBarry Smith const PetscInt *idx; 4847c6ae99SBarry Smith char node[10]; 4947c6ae99SBarry Smith PetscBool isnull; 5047c6ae99SBarry Smith 5147c6ae99SBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 5247c6ae99SBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 536636e97aSMatthew G Knepley if (!da->coordinates) { 5447c6ae99SBarry Smith ierr = PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);CHKERRQ(ierr); 5547c6ae99SBarry Smith } 5647c6ae99SBarry Smith ierr = PetscDrawSynchronizedClear(draw);CHKERRQ(ierr); 5747c6ae99SBarry Smith 5847c6ae99SBarry Smith /* first processor draw all node lines */ 5947c6ae99SBarry Smith if (!rank) { 6047c6ae99SBarry Smith ymin = 0.0; ymax = dd->N - 1; 6147c6ae99SBarry Smith for (xmin=0; xmin<dd->M; xmin++) { 6247c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr); 6347c6ae99SBarry Smith } 6447c6ae99SBarry Smith xmin = 0.0; xmax = dd->M - 1; 6547c6ae99SBarry Smith for (ymin=0; ymin<dd->N; ymin++) { 6647c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);CHKERRQ(ierr); 6747c6ae99SBarry Smith } 6847c6ae99SBarry Smith } 6947c6ae99SBarry Smith ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr); 7047c6ae99SBarry Smith ierr = PetscDrawPause(draw);CHKERRQ(ierr); 7147c6ae99SBarry Smith 7247c6ae99SBarry Smith /* draw my box */ 7347c6ae99SBarry Smith ymin = dd->ys; ymax = dd->ye - 1; xmin = dd->xs/dd->w; 7447c6ae99SBarry Smith xmax =(dd->xe-1)/dd->w; 7547c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);CHKERRQ(ierr); 7647c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 7747c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 7847c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 7947c6ae99SBarry Smith 8047c6ae99SBarry Smith /* put in numbers */ 8147c6ae99SBarry Smith base = (dd->base)/dd->w; 8247c6ae99SBarry Smith for (y=ymin; y<=ymax; y++) { 8347c6ae99SBarry Smith for (x=xmin; x<=xmax; x++) { 8447c6ae99SBarry Smith sprintf(node,"%d",(int)base++); 8547c6ae99SBarry Smith ierr = PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node);CHKERRQ(ierr); 8647c6ae99SBarry Smith } 8747c6ae99SBarry Smith } 8847c6ae99SBarry Smith 8947c6ae99SBarry Smith ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr); 9047c6ae99SBarry Smith ierr = PetscDrawPause(draw);CHKERRQ(ierr); 9147c6ae99SBarry Smith /* overlay ghost numbers, useful for error checking */ 9247c6ae99SBarry Smith /* put in numbers */ 9347c6ae99SBarry Smith 948ea3bf28SBarry Smith base = 0; 9545b6f7e9SBarry Smith ierr = ISLocalToGlobalMappingGetBlockIndices(da->ltogmap,&idx);CHKERRQ(ierr); 9647c6ae99SBarry Smith ymin = dd->Ys; ymax = dd->Ye; xmin = dd->Xs; xmax = dd->Xe; 9747c6ae99SBarry Smith for (y=ymin; y<ymax; y++) { 9847c6ae99SBarry Smith for (x=xmin; x<xmax; x++) { 9947c6ae99SBarry Smith if ((base % dd->w) == 0) { 100bfe97906SBarry Smith sprintf(node,"%d",(int)(idx[base/dd->w])); 10147c6ae99SBarry Smith ierr = PetscDrawString(draw,x/dd->w,y,PETSC_DRAW_BLUE,node);CHKERRQ(ierr); 10247c6ae99SBarry Smith } 10347c6ae99SBarry Smith base++; 10447c6ae99SBarry Smith } 10547c6ae99SBarry Smith } 10645b6f7e9SBarry Smith ierr = ISLocalToGlobalMappingRestoreBlockIndices(da->ltogmap,&idx); 10747c6ae99SBarry Smith ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr); 10847c6ae99SBarry Smith ierr = PetscDrawPause(draw);CHKERRQ(ierr); 1099a42bb27SBarry Smith } else if (isbinary) { 1109a42bb27SBarry Smith ierr = DMView_DA_Binary(da,viewer);CHKERRQ(ierr); 1119a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 1129a42bb27SBarry Smith } else if (ismatlab) { 1139a42bb27SBarry Smith ierr = DMView_DA_Matlab(da,viewer);CHKERRQ(ierr); 1149a42bb27SBarry Smith #endif 11511aeaf0aSBarry Smith } 11647c6ae99SBarry Smith PetscFunctionReturn(0); 11747c6ae99SBarry Smith } 11847c6ae99SBarry Smith 11947c6ae99SBarry Smith /* 12047c6ae99SBarry Smith M is number of grid points 12147c6ae99SBarry Smith m is number of processors 12247c6ae99SBarry Smith 12347c6ae99SBarry Smith */ 12447c6ae99SBarry Smith #undef __FUNCT__ 125aa219208SBarry Smith #define __FUNCT__ "DMDASplitComm2d" 1267087cfbeSBarry Smith PetscErrorCode DMDASplitComm2d(MPI_Comm comm,PetscInt M,PetscInt N,PetscInt sw,MPI_Comm *outcomm) 12747c6ae99SBarry Smith { 12847c6ae99SBarry Smith PetscErrorCode ierr; 12947c6ae99SBarry Smith PetscInt m,n = 0,x = 0,y = 0; 13047c6ae99SBarry Smith PetscMPIInt size,csize,rank; 13147c6ae99SBarry Smith 13247c6ae99SBarry Smith PetscFunctionBegin; 13347c6ae99SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 13447c6ae99SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 13547c6ae99SBarry Smith 13647c6ae99SBarry Smith csize = 4*size; 13747c6ae99SBarry Smith do { 13847c6ae99SBarry 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); 13947c6ae99SBarry Smith csize = csize/4; 14047c6ae99SBarry Smith 141369cc0aeSBarry Smith m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)csize)/((PetscReal)N))); 14247c6ae99SBarry Smith if (!m) m = 1; 14347c6ae99SBarry Smith while (m > 0) { 14447c6ae99SBarry Smith n = csize/m; 14547c6ae99SBarry Smith if (m*n == csize) break; 14647c6ae99SBarry Smith m--; 14747c6ae99SBarry Smith } 14847c6ae99SBarry Smith if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;} 14947c6ae99SBarry Smith 15047c6ae99SBarry Smith x = M/m + ((M % m) > ((csize-1) % m)); 15147c6ae99SBarry Smith y = (N + (csize-1)/m)/n; 15247c6ae99SBarry Smith } while ((x < 4 || y < 4) && csize > 1); 15347c6ae99SBarry Smith if (size != csize) { 15447c6ae99SBarry Smith MPI_Group entire_group,sub_group; 15547c6ae99SBarry Smith PetscMPIInt i,*groupies; 15647c6ae99SBarry Smith 15747c6ae99SBarry Smith ierr = MPI_Comm_group(comm,&entire_group);CHKERRQ(ierr); 158785e854fSJed Brown ierr = PetscMalloc1(csize,&groupies);CHKERRQ(ierr); 15947c6ae99SBarry Smith for (i=0; i<csize; i++) { 16047c6ae99SBarry Smith groupies[i] = (rank/csize)*csize + i; 16147c6ae99SBarry Smith } 16247c6ae99SBarry Smith ierr = MPI_Group_incl(entire_group,csize,groupies,&sub_group);CHKERRQ(ierr); 16347c6ae99SBarry Smith ierr = PetscFree(groupies);CHKERRQ(ierr); 16447c6ae99SBarry Smith ierr = MPI_Comm_create(comm,sub_group,outcomm);CHKERRQ(ierr); 16547c6ae99SBarry Smith ierr = MPI_Group_free(&entire_group);CHKERRQ(ierr); 16647c6ae99SBarry Smith ierr = MPI_Group_free(&sub_group);CHKERRQ(ierr); 167aa219208SBarry Smith ierr = PetscInfo1(0,"DMDASplitComm2d:Creating redundant coarse problems of size %d\n",csize);CHKERRQ(ierr); 16847c6ae99SBarry Smith } else { 16947c6ae99SBarry Smith *outcomm = comm; 17047c6ae99SBarry Smith } 17147c6ae99SBarry Smith PetscFunctionReturn(0); 17247c6ae99SBarry Smith } 17347c6ae99SBarry Smith 17447c6ae99SBarry Smith #if defined(new) 17547c6ae99SBarry Smith #undef __FUNCT__ 176aa219208SBarry Smith #define __FUNCT__ "DMDAGetDiagonal_MFFD" 17747c6ae99SBarry Smith /* 178aa219208SBarry Smith DMDAGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix where local 179aa219208SBarry Smith function lives on a DMDA 18047c6ae99SBarry Smith 18147c6ae99SBarry Smith y ~= (F(u + ha) - F(u))/h, 18247c6ae99SBarry Smith where F = nonlinear function, as set by SNESSetFunction() 18347c6ae99SBarry Smith u = current iterate 18447c6ae99SBarry Smith h = difference interval 18547c6ae99SBarry Smith */ 186aa219208SBarry Smith PetscErrorCode DMDAGetDiagonal_MFFD(DM da,Vec U,Vec a) 18747c6ae99SBarry Smith { 18847c6ae99SBarry Smith PetscScalar h,*aa,*ww,v; 18947c6ae99SBarry Smith PetscReal epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON; 19047c6ae99SBarry Smith PetscErrorCode ierr; 19147c6ae99SBarry Smith PetscInt gI,nI; 19247c6ae99SBarry Smith MatStencil stencil; 193aa219208SBarry Smith DMDALocalInfo info; 19447c6ae99SBarry Smith 19547c6ae99SBarry Smith PetscFunctionBegin; 19647c6ae99SBarry Smith ierr = (*ctx->func)(0,U,a,ctx->funcctx);CHKERRQ(ierr); 19747c6ae99SBarry Smith ierr = (*ctx->funcisetbase)(U,ctx->funcctx);CHKERRQ(ierr); 19847c6ae99SBarry Smith 19947c6ae99SBarry Smith ierr = VecGetArray(U,&ww);CHKERRQ(ierr); 20047c6ae99SBarry Smith ierr = VecGetArray(a,&aa);CHKERRQ(ierr); 20147c6ae99SBarry Smith 20247c6ae99SBarry Smith nI = 0; 20347c6ae99SBarry Smith h = ww[gI]; 20447c6ae99SBarry Smith if (h == 0.0) h = 1.0; 20547c6ae99SBarry Smith if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin; 20647c6ae99SBarry Smith else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin; 20747c6ae99SBarry Smith h *= epsilon; 20847c6ae99SBarry Smith 20947c6ae99SBarry Smith ww[gI] += h; 21047c6ae99SBarry Smith ierr = (*ctx->funci)(i,w,&v,ctx->funcctx);CHKERRQ(ierr); 21147c6ae99SBarry Smith aa[nI] = (v - aa[nI])/h; 21247c6ae99SBarry Smith ww[gI] -= h; 21347c6ae99SBarry Smith nI++; 2148865f1eaSKarl Rupp 21547c6ae99SBarry Smith ierr = VecRestoreArray(U,&ww);CHKERRQ(ierr); 21647c6ae99SBarry Smith ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr); 21747c6ae99SBarry Smith PetscFunctionReturn(0); 21847c6ae99SBarry Smith } 21947c6ae99SBarry Smith #endif 22047c6ae99SBarry Smith 22147c6ae99SBarry Smith #undef __FUNCT__ 2229a42bb27SBarry Smith #define __FUNCT__ "DMSetUp_DA_2D" 2237087cfbeSBarry Smith PetscErrorCode DMSetUp_DA_2D(DM da) 22447c6ae99SBarry Smith { 22547c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 22647c6ae99SBarry Smith const PetscInt M = dd->M; 22747c6ae99SBarry Smith const PetscInt N = dd->N; 22847c6ae99SBarry Smith PetscInt m = dd->m; 22947c6ae99SBarry Smith PetscInt n = dd->n; 23047c6ae99SBarry Smith const PetscInt dof = dd->w; 23147c6ae99SBarry Smith const PetscInt s = dd->s; 232bff4a2f0SMatthew G. Knepley DMBoundaryType bx = dd->bx; 233bff4a2f0SMatthew G. Knepley DMBoundaryType by = dd->by; 23419fd82e9SBarry Smith DMDAStencilType stencil_type = dd->stencil_type; 23547c6ae99SBarry Smith PetscInt *lx = dd->lx; 23647c6ae99SBarry Smith PetscInt *ly = dd->ly; 23747c6ae99SBarry Smith MPI_Comm comm; 23847c6ae99SBarry Smith PetscMPIInt rank,size; 239ce00eea3SSatish Balay PetscInt xs,xe,ys,ye,x,y,Xs,Xe,Ys,Ye,start,end,IXs,IXe,IYs,IYe; 2408ea3bf28SBarry Smith PetscInt up,down,left,right,i,n0,n1,n2,n3,n5,n6,n7,n8,*idx,nn; 241db87c5ecSEthan Coon PetscInt xbase,*bases,*ldims,j,x_t,y_t,s_t,base,count; 24247c6ae99SBarry Smith PetscInt s_x,s_y; /* s proportionalized to w */ 24347c6ae99SBarry Smith PetscInt sn0 = 0,sn2 = 0,sn6 = 0,sn8 = 0; 24447c6ae99SBarry Smith Vec local,global; 24547c6ae99SBarry Smith VecScatter ltog,gtol; 24645b6f7e9SBarry Smith IS to,from; 24747c6ae99SBarry Smith PetscErrorCode ierr; 24847c6ae99SBarry Smith 24947c6ae99SBarry Smith PetscFunctionBegin; 250bff4a2f0SMatthew 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"); 25147c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 2523855c12bSBarry Smith #if !defined(PETSC_USE_64BIT_INDICES) 2533855c12bSBarry Smith if (((Petsc64bitInt) M)*((Petsc64bitInt) N)*((Petsc64bitInt) dof) > (Petsc64bitInt) 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); 2543855c12bSBarry Smith #endif 25547c6ae99SBarry Smith 25647c6ae99SBarry Smith if (dof < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have 1 or more degrees of freedom per node: %D",dof); 25747c6ae99SBarry Smith if (s < 0) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Stencil width cannot be negative: %D",s); 25847c6ae99SBarry Smith 25947c6ae99SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 26047c6ae99SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 26147c6ae99SBarry Smith 26247c6ae99SBarry Smith if (m != PETSC_DECIDE) { 26347c6ae99SBarry Smith if (m < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m); 26447c6ae99SBarry Smith else if (m > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size); 26547c6ae99SBarry Smith } 26647c6ae99SBarry Smith if (n != PETSC_DECIDE) { 26747c6ae99SBarry Smith if (n < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n); 26847c6ae99SBarry Smith else if (n > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size); 26947c6ae99SBarry Smith } 27047c6ae99SBarry Smith 27147c6ae99SBarry Smith if (m == PETSC_DECIDE || n == PETSC_DECIDE) { 27247c6ae99SBarry Smith if (n != PETSC_DECIDE) { 27347c6ae99SBarry Smith m = size/n; 27447c6ae99SBarry Smith } else if (m != PETSC_DECIDE) { 27547c6ae99SBarry Smith n = size/m; 27647c6ae99SBarry Smith } else { 27747c6ae99SBarry Smith /* try for squarish distribution */ 278369cc0aeSBarry Smith m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)N))); 27947c6ae99SBarry Smith if (!m) m = 1; 28047c6ae99SBarry Smith while (m > 0) { 28147c6ae99SBarry Smith n = size/m; 28247c6ae99SBarry Smith if (m*n == size) break; 28347c6ae99SBarry Smith m--; 28447c6ae99SBarry Smith } 28547c6ae99SBarry Smith if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;} 28647c6ae99SBarry Smith } 28747c6ae99SBarry 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 "); 28847c6ae99SBarry Smith } else if (m*n != size) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition"); 28947c6ae99SBarry Smith 29047c6ae99SBarry Smith if (M < m) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m); 29147c6ae99SBarry Smith if (N < n) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n); 29247c6ae99SBarry Smith 29347c6ae99SBarry Smith /* 29447c6ae99SBarry Smith Determine locally owned region 29547c6ae99SBarry Smith xs is the first local node number, x is the number of local nodes 29647c6ae99SBarry Smith */ 29747c6ae99SBarry Smith if (!lx) { 298785e854fSJed Brown ierr = PetscMalloc1(m, &dd->lx);CHKERRQ(ierr); 29947c6ae99SBarry Smith lx = dd->lx; 30047c6ae99SBarry Smith for (i=0; i<m; i++) { 30147c6ae99SBarry Smith lx[i] = M/m + ((M % m) > i); 30247c6ae99SBarry Smith } 30347c6ae99SBarry Smith } 30447c6ae99SBarry Smith x = lx[rank % m]; 30547c6ae99SBarry Smith xs = 0; 30647c6ae99SBarry Smith for (i=0; i<(rank % m); i++) { 30747c6ae99SBarry Smith xs += lx[i]; 30847c6ae99SBarry Smith } 30947c6ae99SBarry Smith #if defined(PETSC_USE_DEBUG) 31047c6ae99SBarry Smith left = xs; 31147c6ae99SBarry Smith for (i=(rank % m); i<m; i++) { 31247c6ae99SBarry Smith left += lx[i]; 31347c6ae99SBarry Smith } 31447c6ae99SBarry 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); 31547c6ae99SBarry Smith #endif 31647c6ae99SBarry Smith 31747c6ae99SBarry Smith /* 31847c6ae99SBarry Smith Determine locally owned region 31947c6ae99SBarry Smith ys is the first local node number, y is the number of local nodes 32047c6ae99SBarry Smith */ 32147c6ae99SBarry Smith if (!ly) { 322785e854fSJed Brown ierr = PetscMalloc1(n, &dd->ly);CHKERRQ(ierr); 32347c6ae99SBarry Smith ly = dd->ly; 32447c6ae99SBarry Smith for (i=0; i<n; i++) { 32547c6ae99SBarry Smith ly[i] = N/n + ((N % n) > i); 32647c6ae99SBarry Smith } 32747c6ae99SBarry Smith } 32847c6ae99SBarry Smith y = ly[rank/m]; 32947c6ae99SBarry Smith ys = 0; 33047c6ae99SBarry Smith for (i=0; i<(rank/m); i++) { 33147c6ae99SBarry Smith ys += ly[i]; 33247c6ae99SBarry Smith } 33347c6ae99SBarry Smith #if defined(PETSC_USE_DEBUG) 33447c6ae99SBarry Smith left = ys; 33547c6ae99SBarry Smith for (i=(rank/m); i<n; i++) { 33647c6ae99SBarry Smith left += ly[i]; 33747c6ae99SBarry Smith } 33847c6ae99SBarry 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); 33947c6ae99SBarry Smith #endif 34047c6ae99SBarry Smith 341bcea557cSEthan Coon /* 342bcea557cSEthan Coon check if the scatter requires more than one process neighbor or wraps around 343bcea557cSEthan Coon the domain more than once 344bcea557cSEthan Coon */ 345bff4a2f0SMatthew 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); 346bff4a2f0SMatthew 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); 34747c6ae99SBarry Smith xe = xs + x; 34847c6ae99SBarry Smith ye = ys + y; 34947c6ae99SBarry Smith 350ce00eea3SSatish Balay /* determine ghost region (Xs) and region scattered into (IXs) */ 351d9c9ebe5SPeter Brune if (xs-s > 0) { 352d9c9ebe5SPeter Brune Xs = xs - s; IXs = xs - s; 35388661749SPeter Brune } else { 35488661749SPeter Brune if (bx) { 35588661749SPeter Brune Xs = xs - s; 35688661749SPeter Brune } else { 35788661749SPeter Brune Xs = 0; 35888661749SPeter Brune } 35988661749SPeter Brune IXs = 0; 36088661749SPeter Brune } 361d9c9ebe5SPeter Brune if (xe+s <= M) { 362d9c9ebe5SPeter Brune Xe = xe + s; IXe = xe + s; 36388661749SPeter Brune } else { 36488661749SPeter Brune if (bx) { 365d9c9ebe5SPeter Brune Xs = xs - s; Xe = xe + s; 36688661749SPeter Brune } else { 36788661749SPeter Brune Xe = M; 36888661749SPeter Brune } 36988661749SPeter Brune IXe = M; 37088661749SPeter Brune } 37147c6ae99SBarry Smith 372bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC || bx == DM_BOUNDARY_MIRROR) { 373d9c9ebe5SPeter Brune IXs = xs - s; 374d9c9ebe5SPeter Brune IXe = xe + s; 375d9c9ebe5SPeter Brune Xs = xs - s; 376d9c9ebe5SPeter Brune Xe = xe + s; 37788661749SPeter Brune } 37847c6ae99SBarry Smith 379d9c9ebe5SPeter Brune if (ys-s > 0) { 380d9c9ebe5SPeter Brune Ys = ys - s; IYs = ys - s; 38188661749SPeter Brune } else { 38288661749SPeter Brune if (by) { 38388661749SPeter Brune Ys = ys - s; 38488661749SPeter Brune } else { 38588661749SPeter Brune Ys = 0; 38688661749SPeter Brune } 38788661749SPeter Brune IYs = 0; 38888661749SPeter Brune } 389d9c9ebe5SPeter Brune if (ye+s <= N) { 390d9c9ebe5SPeter Brune Ye = ye + s; IYe = ye + s; 39188661749SPeter Brune } else { 39288661749SPeter Brune if (by) { 39388661749SPeter Brune Ye = ye + s; 39488661749SPeter Brune } else { 39588661749SPeter Brune Ye = N; 39688661749SPeter Brune } 39788661749SPeter Brune IYe = N; 39888661749SPeter Brune } 39988661749SPeter Brune 400bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_PERIODIC || by == DM_BOUNDARY_MIRROR) { 401d9c9ebe5SPeter Brune IYs = ys - s; 402d9c9ebe5SPeter Brune IYe = ye + s; 403d9c9ebe5SPeter Brune Ys = ys - s; 404d9c9ebe5SPeter Brune Ye = ye + s; 40588661749SPeter Brune } 40688661749SPeter Brune 40788661749SPeter Brune /* stencil length in each direction */ 408d9c9ebe5SPeter Brune s_x = s; 409d9c9ebe5SPeter Brune s_y = s; 41047c6ae99SBarry Smith 41147c6ae99SBarry Smith /* determine starting point of each processor */ 41247c6ae99SBarry Smith nn = x*y; 413dcca6d9dSJed Brown ierr = PetscMalloc2(size+1,&bases,size,&ldims);CHKERRQ(ierr); 41447c6ae99SBarry Smith ierr = MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);CHKERRQ(ierr); 41547c6ae99SBarry Smith bases[0] = 0; 41647c6ae99SBarry Smith for (i=1; i<=size; i++) { 41747c6ae99SBarry Smith bases[i] = ldims[i-1]; 41847c6ae99SBarry Smith } 41947c6ae99SBarry Smith for (i=1; i<=size; i++) { 42047c6ae99SBarry Smith bases[i] += bases[i-1]; 42147c6ae99SBarry Smith } 422ce00eea3SSatish Balay base = bases[rank]*dof; 42347c6ae99SBarry Smith 42447c6ae99SBarry Smith /* allocate the base parallel and sequential vectors */ 425ce00eea3SSatish Balay dd->Nlocal = x*y*dof; 426b1fb7eb7SBarry Smith ierr = VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,NULL,&global);CHKERRQ(ierr); 427ce00eea3SSatish Balay dd->nlocal = (Xe-Xs)*(Ye-Ys)*dof; 428b1fb7eb7SBarry Smith ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,NULL,&local);CHKERRQ(ierr); 42947c6ae99SBarry Smith 43047c6ae99SBarry Smith /* generate appropriate vector scatters */ 43147c6ae99SBarry Smith /* local to global inserts non-ghost point region into global */ 43247c6ae99SBarry Smith ierr = VecGetOwnershipRange(global,&start,&end);CHKERRQ(ierr); 433ce00eea3SSatish Balay ierr = ISCreateStride(comm,x*y*dof,start,1,&to);CHKERRQ(ierr); 43447c6ae99SBarry Smith 435*02fe608eSBarry Smith ierr = PetscMalloc1((IXe-IXs)*(IYe-IYs),&idx);CHKERRQ(ierr); 436ce00eea3SSatish Balay left = xs - Xs; right = left + x; 437ce00eea3SSatish Balay down = ys - Ys; up = down + y; 43847c6ae99SBarry Smith count = 0; 43947c6ae99SBarry Smith for (i=down; i<up; i++) { 440ce00eea3SSatish Balay for (j=left; j<right; j++) { 441ce00eea3SSatish Balay idx[count++] = i*(Xe-Xs) + j; 44247c6ae99SBarry Smith } 44347c6ae99SBarry Smith } 44447c6ae99SBarry Smith 445*02fe608eSBarry Smith ierr = ISCreateBlock(comm,dof,count,idx,PETSC_USE_POINTER,&from);CHKERRQ(ierr); 44647c6ae99SBarry Smith ierr = VecScatterCreate(local,from,global,to,<og);CHKERRQ(ierr); 447f4a8f317SJed Brown ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)ltog);CHKERRQ(ierr); 448fcfd50ebSBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 449fcfd50ebSBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 45047c6ae99SBarry Smith 451ce00eea3SSatish Balay /* global to local must include ghost points within the domain, 452ce00eea3SSatish Balay but not ghost points outside the domain that aren't periodic */ 453d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_BOX) { 454ce00eea3SSatish Balay left = IXs - Xs; right = left + (IXe-IXs); 455ce00eea3SSatish Balay down = IYs - Ys; up = down + (IYe-IYs); 456ce00eea3SSatish Balay count = 0; 457ce00eea3SSatish Balay for (i=down; i<up; i++) { 458ce00eea3SSatish Balay for (j=left; j<right; j++) { 459ce00eea3SSatish Balay idx[count++] = j + i*(Xe-Xs); 460ce00eea3SSatish Balay } 461ce00eea3SSatish Balay } 462ce00eea3SSatish Balay ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr); 463ce00eea3SSatish Balay 46447c6ae99SBarry Smith } else { 46547c6ae99SBarry Smith /* must drop into cross shape region */ 46647c6ae99SBarry Smith /* ---------| 46747c6ae99SBarry Smith | top | 468ce00eea3SSatish Balay |--- ---| up 46947c6ae99SBarry Smith | middle | 47047c6ae99SBarry Smith | | 471ce00eea3SSatish Balay ---- ---- down 47247c6ae99SBarry Smith | bottom | 47347c6ae99SBarry Smith ----------- 47447c6ae99SBarry Smith Xs xs xe Xe */ 475ce00eea3SSatish Balay left = xs - Xs; right = left + x; 476ce00eea3SSatish Balay down = ys - Ys; up = down + y; 47747c6ae99SBarry Smith count = 0; 478ce00eea3SSatish Balay /* bottom */ 479ce00eea3SSatish Balay for (i=(IYs-Ys); i<down; i++) { 480ce00eea3SSatish Balay for (j=left; j<right; j++) { 481ce00eea3SSatish Balay idx[count++] = j + i*(Xe-Xs); 48247c6ae99SBarry Smith } 48347c6ae99SBarry Smith } 48447c6ae99SBarry Smith /* middle */ 48547c6ae99SBarry Smith for (i=down; i<up; i++) { 486ce00eea3SSatish Balay for (j=(IXs-Xs); j<(IXe-Xs); j++) { 487ce00eea3SSatish Balay idx[count++] = j + i*(Xe-Xs); 48847c6ae99SBarry Smith } 48947c6ae99SBarry Smith } 49047c6ae99SBarry Smith /* top */ 491ce00eea3SSatish Balay for (i=up; i<up+IYe-ye; i++) { 492ce00eea3SSatish Balay for (j=left; j<right; j++) { 493ce00eea3SSatish Balay idx[count++] = j + i*(Xe-Xs); 49447c6ae99SBarry Smith } 49547c6ae99SBarry Smith } 49647c6ae99SBarry Smith ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr); 49747c6ae99SBarry Smith } 49847c6ae99SBarry Smith 49947c6ae99SBarry Smith 50047c6ae99SBarry Smith /* determine who lies on each side of us stored in n6 n7 n8 50147c6ae99SBarry Smith n3 n5 50247c6ae99SBarry Smith n0 n1 n2 50347c6ae99SBarry Smith */ 50447c6ae99SBarry Smith 50547c6ae99SBarry Smith /* Assume the Non-Periodic Case */ 50647c6ae99SBarry Smith n1 = rank - m; 50747c6ae99SBarry Smith if (rank % m) { 50847c6ae99SBarry Smith n0 = n1 - 1; 50947c6ae99SBarry Smith } else { 51047c6ae99SBarry Smith n0 = -1; 51147c6ae99SBarry Smith } 51247c6ae99SBarry Smith if ((rank+1) % m) { 51347c6ae99SBarry Smith n2 = n1 + 1; 51447c6ae99SBarry Smith n5 = rank + 1; 51547c6ae99SBarry Smith n8 = rank + m + 1; if (n8 >= m*n) n8 = -1; 51647c6ae99SBarry Smith } else { 51747c6ae99SBarry Smith n2 = -1; n5 = -1; n8 = -1; 51847c6ae99SBarry Smith } 51947c6ae99SBarry Smith if (rank % m) { 52047c6ae99SBarry Smith n3 = rank - 1; 52147c6ae99SBarry Smith n6 = n3 + m; if (n6 >= m*n) n6 = -1; 52247c6ae99SBarry Smith } else { 52347c6ae99SBarry Smith n3 = -1; n6 = -1; 52447c6ae99SBarry Smith } 52547c6ae99SBarry Smith n7 = rank + m; if (n7 >= m*n) n7 = -1; 52647c6ae99SBarry Smith 527bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC && by == DM_BOUNDARY_PERIODIC) { 52847c6ae99SBarry Smith /* Modify for Periodic Cases */ 52947c6ae99SBarry Smith /* Handle all four corners */ 53047c6ae99SBarry Smith if ((n6 < 0) && (n7 < 0) && (n3 < 0)) n6 = m-1; 53147c6ae99SBarry Smith if ((n8 < 0) && (n7 < 0) && (n5 < 0)) n8 = 0; 53247c6ae99SBarry Smith if ((n2 < 0) && (n5 < 0) && (n1 < 0)) n2 = size-m; 53347c6ae99SBarry Smith if ((n0 < 0) && (n3 < 0) && (n1 < 0)) n0 = size-1; 53447c6ae99SBarry Smith 53547c6ae99SBarry Smith /* Handle Top and Bottom Sides */ 53647c6ae99SBarry Smith if (n1 < 0) n1 = rank + m * (n-1); 53747c6ae99SBarry Smith if (n7 < 0) n7 = rank - m * (n-1); 53847c6ae99SBarry Smith if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1; 53947c6ae99SBarry Smith if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1; 54047c6ae99SBarry Smith if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1; 54147c6ae99SBarry Smith if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1; 54247c6ae99SBarry Smith 54347c6ae99SBarry Smith /* Handle Left and Right Sides */ 54447c6ae99SBarry Smith if (n3 < 0) n3 = rank + (m-1); 54547c6ae99SBarry Smith if (n5 < 0) n5 = rank - (m-1); 54647c6ae99SBarry Smith if ((n1 >= 0) && (n0 < 0)) n0 = rank-1; 54747c6ae99SBarry Smith if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1; 54847c6ae99SBarry Smith if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1; 54947c6ae99SBarry Smith if ((n7 >= 0) && (n8 < 0)) n8 = rank+1; 550bff4a2f0SMatthew G. Knepley } else if (by == DM_BOUNDARY_PERIODIC) { /* Handle Top and Bottom Sides */ 551ce00eea3SSatish Balay if (n1 < 0) n1 = rank + m * (n-1); 552ce00eea3SSatish Balay if (n7 < 0) n7 = rank - m * (n-1); 553ce00eea3SSatish Balay if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1; 554ce00eea3SSatish Balay if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1; 555ce00eea3SSatish Balay if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1; 556ce00eea3SSatish Balay if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1; 557bff4a2f0SMatthew G. Knepley } else if (bx == DM_BOUNDARY_PERIODIC) { /* Handle Left and Right Sides */ 558ce00eea3SSatish Balay if (n3 < 0) n3 = rank + (m-1); 559ce00eea3SSatish Balay if (n5 < 0) n5 = rank - (m-1); 560ce00eea3SSatish Balay if ((n1 >= 0) && (n0 < 0)) n0 = rank-1; 561ce00eea3SSatish Balay if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1; 562ce00eea3SSatish Balay if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1; 563ce00eea3SSatish Balay if ((n7 >= 0) && (n8 < 0)) n8 = rank+1; 56447c6ae99SBarry Smith } 565ce00eea3SSatish Balay 566785e854fSJed Brown ierr = PetscMalloc1(9,&dd->neighbors);CHKERRQ(ierr); 5678865f1eaSKarl Rupp 56847c6ae99SBarry Smith dd->neighbors[0] = n0; 56947c6ae99SBarry Smith dd->neighbors[1] = n1; 57047c6ae99SBarry Smith dd->neighbors[2] = n2; 57147c6ae99SBarry Smith dd->neighbors[3] = n3; 57247c6ae99SBarry Smith dd->neighbors[4] = rank; 57347c6ae99SBarry Smith dd->neighbors[5] = n5; 57447c6ae99SBarry Smith dd->neighbors[6] = n6; 57547c6ae99SBarry Smith dd->neighbors[7] = n7; 57647c6ae99SBarry Smith dd->neighbors[8] = n8; 57747c6ae99SBarry Smith 578d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_STAR) { 57947c6ae99SBarry Smith /* save corner processor numbers */ 58047c6ae99SBarry Smith sn0 = n0; sn2 = n2; sn6 = n6; sn8 = n8; 58147c6ae99SBarry Smith n0 = n2 = n6 = n8 = -1; 58247c6ae99SBarry Smith } 58347c6ae99SBarry Smith 584785e854fSJed Brown ierr = PetscMalloc1((Xe-Xs)*(Ye-Ys),&idx);CHKERRQ(ierr); 5853bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)da,(Xe-Xs)*(Ye-Ys)*sizeof(PetscInt));CHKERRQ(ierr); 58647c6ae99SBarry Smith 587ce00eea3SSatish Balay nn = 0; 58847c6ae99SBarry Smith xbase = bases[rank]; 58947c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 59047c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 591ce00eea3SSatish Balay x_t = lx[n0 % m]; 59247c6ae99SBarry Smith y_t = ly[(n0/m)]; 59347c6ae99SBarry Smith s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x; 5948865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 59547c6ae99SBarry Smith } 596ac119b13SBarry Smith 59747c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 59847c6ae99SBarry Smith x_t = x; 59947c6ae99SBarry Smith y_t = ly[(n1/m)]; 60047c6ae99SBarry Smith s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t; 6018865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 602bff4a2f0SMatthew G. Knepley } else if (by == DM_BOUNDARY_MIRROR) { 6038865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(s_y - i + 1) + j; 60447c6ae99SBarry Smith } 605ac119b13SBarry Smith 60647c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 607ce00eea3SSatish Balay x_t = lx[n2 % m]; 60847c6ae99SBarry Smith y_t = ly[(n2/m)]; 60947c6ae99SBarry Smith s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t; 6108865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 61147c6ae99SBarry Smith } 61247c6ae99SBarry Smith } 61347c6ae99SBarry Smith 61447c6ae99SBarry Smith for (i=0; i<y; i++) { 61547c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 616ce00eea3SSatish Balay x_t = lx[n3 % m]; 61747c6ae99SBarry Smith /* y_t = y; */ 61847c6ae99SBarry Smith s_t = bases[n3] + (i+1)*x_t - s_x; 6198865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 620bff4a2f0SMatthew G. Knepley } else if (bx == DM_BOUNDARY_MIRROR) { 6218865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*i + s_x - j; 62247c6ae99SBarry Smith } 62347c6ae99SBarry Smith 6248865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = xbase++; /* interior */ 62547c6ae99SBarry Smith 62647c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 627ce00eea3SSatish Balay x_t = lx[n5 % m]; 62847c6ae99SBarry Smith /* y_t = y; */ 62947c6ae99SBarry Smith s_t = bases[n5] + (i)*x_t; 6308865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 631bff4a2f0SMatthew G. Knepley } else if (bx == DM_BOUNDARY_MIRROR) { 6328865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*(i + 1) - 2 - j; 63347c6ae99SBarry Smith } 63447c6ae99SBarry Smith } 63547c6ae99SBarry Smith 63647c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 63747c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 638ce00eea3SSatish Balay x_t = lx[n6 % m]; 63947c6ae99SBarry Smith /* y_t = ly[(n6/m)]; */ 64047c6ae99SBarry Smith s_t = bases[n6] + (i)*x_t - s_x; 6418865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 64247c6ae99SBarry Smith } 643ac119b13SBarry Smith 64447c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 64547c6ae99SBarry Smith x_t = x; 64647c6ae99SBarry Smith /* y_t = ly[(n7/m)]; */ 64747c6ae99SBarry Smith s_t = bases[n7] + (i-1)*x_t; 6488865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 649bff4a2f0SMatthew G. Knepley } else if (by == DM_BOUNDARY_MIRROR) { 6508865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(y - i - 1) + j; 65147c6ae99SBarry Smith } 652ac119b13SBarry Smith 65347c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 654ce00eea3SSatish Balay x_t = lx[n8 % m]; 65547c6ae99SBarry Smith /* y_t = ly[(n8/m)]; */ 65647c6ae99SBarry Smith s_t = bases[n8] + (i-1)*x_t; 6578865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 65847c6ae99SBarry Smith } 65947c6ae99SBarry Smith } 66047c6ae99SBarry Smith 661b1fb7eb7SBarry Smith ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_USE_POINTER,&from);CHKERRQ(ierr); 66247c6ae99SBarry Smith ierr = VecScatterCreate(global,from,local,to,>ol);CHKERRQ(ierr); 6633bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)gtol);CHKERRQ(ierr); 664fcfd50ebSBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 665fcfd50ebSBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 66647c6ae99SBarry Smith 667d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_STAR) { 668ce00eea3SSatish Balay n0 = sn0; n2 = sn2; n6 = sn6; n8 = sn8; 669ce00eea3SSatish Balay } 670ce00eea3SSatish Balay 67188661749SPeter Brune if (((stencil_type == DMDA_STENCIL_STAR) || 672bff4a2f0SMatthew G. Knepley (bx && bx != DM_BOUNDARY_PERIODIC) || 673bff4a2f0SMatthew G. Knepley (by && by != DM_BOUNDARY_PERIODIC))) { 67447c6ae99SBarry Smith /* 67547c6ae99SBarry Smith Recompute the local to global mappings, this time keeping the 676ce00eea3SSatish Balay information about the cross corner processor numbers and any ghosted 677ce00eea3SSatish Balay but not periodic indices. 67847c6ae99SBarry Smith */ 67947c6ae99SBarry Smith nn = 0; 68047c6ae99SBarry Smith xbase = bases[rank]; 68147c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 68247c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 683ce00eea3SSatish Balay x_t = lx[n0 % m]; 68447c6ae99SBarry Smith y_t = ly[(n0/m)]; 68547c6ae99SBarry Smith s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x; 6868865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 687ce00eea3SSatish Balay } else if (xs-Xs > 0 && ys-Ys > 0) { 6888865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 68947c6ae99SBarry Smith } 69047c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 69147c6ae99SBarry Smith x_t = x; 69247c6ae99SBarry Smith y_t = ly[(n1/m)]; 69347c6ae99SBarry Smith s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t; 6948865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 695ce00eea3SSatish Balay } else if (ys-Ys > 0) { 696bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_MIRROR) { 6978865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(s_y - i + 1) + j; 698624904c4SBarry Smith } else { 6998865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 70047c6ae99SBarry Smith } 701624904c4SBarry Smith } 70247c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 703ce00eea3SSatish Balay x_t = lx[n2 % m]; 70447c6ae99SBarry Smith y_t = ly[(n2/m)]; 70547c6ae99SBarry Smith s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t; 7068865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 707ce00eea3SSatish Balay } else if (Xe-xe> 0 && ys-Ys > 0) { 7088865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 70947c6ae99SBarry Smith } 71047c6ae99SBarry Smith } 71147c6ae99SBarry Smith 71247c6ae99SBarry Smith for (i=0; i<y; i++) { 71347c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 714ce00eea3SSatish Balay x_t = lx[n3 % m]; 71547c6ae99SBarry Smith /* y_t = y; */ 71647c6ae99SBarry Smith s_t = bases[n3] + (i+1)*x_t - s_x; 7178865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 718ce00eea3SSatish Balay } else if (xs-Xs > 0) { 719bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_MIRROR) { 7208865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*i + s_x - j; 721624904c4SBarry Smith } else { 7228865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 72347c6ae99SBarry Smith } 724624904c4SBarry Smith } 72547c6ae99SBarry Smith 7268865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = xbase++; /* interior */ 72747c6ae99SBarry Smith 72847c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 729ce00eea3SSatish Balay x_t = lx[n5 % m]; 73047c6ae99SBarry Smith /* y_t = y; */ 73147c6ae99SBarry Smith s_t = bases[n5] + (i)*x_t; 7328865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 733ce00eea3SSatish Balay } else if (Xe-xe > 0) { 734bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_MIRROR) { 7358865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*(i + 1) - 2 - j; 736624904c4SBarry Smith } else { 7378865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 73847c6ae99SBarry Smith } 73947c6ae99SBarry Smith } 740624904c4SBarry Smith } 74147c6ae99SBarry Smith 74247c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 74347c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 744ce00eea3SSatish Balay x_t = lx[n6 % m]; 74547c6ae99SBarry Smith /* y_t = ly[(n6/m)]; */ 74647c6ae99SBarry Smith s_t = bases[n6] + (i)*x_t - s_x; 7478865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 748ce00eea3SSatish Balay } else if (xs-Xs > 0 && Ye-ye > 0) { 7498865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 75047c6ae99SBarry Smith } 75147c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 75247c6ae99SBarry Smith x_t = x; 75347c6ae99SBarry Smith /* y_t = ly[(n7/m)]; */ 75447c6ae99SBarry Smith s_t = bases[n7] + (i-1)*x_t; 7558865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 756ce00eea3SSatish Balay } else if (Ye-ye > 0) { 757bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_MIRROR) { 7588865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(y - i - 1) + j; 759624904c4SBarry Smith } else { 7608865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 76147c6ae99SBarry Smith } 762624904c4SBarry Smith } 76347c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 764ce00eea3SSatish Balay x_t = lx[n8 % m]; 76547c6ae99SBarry Smith /* y_t = ly[(n8/m)]; */ 76647c6ae99SBarry Smith s_t = bases[n8] + (i-1)*x_t; 7678865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 768ce00eea3SSatish Balay } else if (Xe-xe > 0 && Ye-ye > 0) { 7698865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 77047c6ae99SBarry Smith } 77147c6ae99SBarry Smith } 77247c6ae99SBarry Smith } 773ce00eea3SSatish Balay /* 774ce00eea3SSatish Balay Set the local to global ordering in the global vector, this allows use 775ce00eea3SSatish Balay of VecSetValuesLocal(). 776ce00eea3SSatish Balay */ 77745b6f7e9SBarry Smith ierr = ISLocalToGlobalMappingCreate(comm,dof,nn,idx,PETSC_OWN_POINTER,&da->ltogmap);CHKERRQ(ierr); 7783bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)da->ltogmap);CHKERRQ(ierr); 77947c6ae99SBarry Smith 780ce00eea3SSatish Balay ierr = PetscFree2(bases,ldims);CHKERRQ(ierr); 78147c6ae99SBarry Smith dd->m = m; dd->n = n; 782ce00eea3SSatish Balay /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */ 783ce00eea3SSatish Balay dd->xs = xs*dof; dd->xe = xe*dof; dd->ys = ys; dd->ye = ye; dd->zs = 0; dd->ze = 1; 784ce00eea3SSatish Balay dd->Xs = Xs*dof; dd->Xe = Xe*dof; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = 0; dd->Ze = 1; 78547c6ae99SBarry Smith 786fcfd50ebSBarry Smith ierr = VecDestroy(&local);CHKERRQ(ierr); 787fcfd50ebSBarry Smith ierr = VecDestroy(&global);CHKERRQ(ierr); 78847c6ae99SBarry Smith 78947c6ae99SBarry Smith dd->gtol = gtol; 79047c6ae99SBarry Smith dd->ltog = ltog; 79147c6ae99SBarry Smith dd->base = base; 7929a42bb27SBarry Smith da->ops->view = DMView_DA_2d; 7930298fd71SBarry Smith dd->ltol = NULL; 7940298fd71SBarry Smith dd->ao = NULL; 79547c6ae99SBarry Smith PetscFunctionReturn(0); 79647c6ae99SBarry Smith } 79747c6ae99SBarry Smith 79847c6ae99SBarry Smith #undef __FUNCT__ 799aa219208SBarry Smith #define __FUNCT__ "DMDACreate2d" 80047c6ae99SBarry Smith /*@C 801aa219208SBarry Smith DMDACreate2d - Creates an object that will manage the communication of two-dimensional 80247c6ae99SBarry Smith regular array data that is distributed across some processors. 80347c6ae99SBarry Smith 80447c6ae99SBarry Smith Collective on MPI_Comm 80547c6ae99SBarry Smith 80647c6ae99SBarry Smith Input Parameters: 80747c6ae99SBarry Smith + comm - MPI communicator 8081321219cSEthan Coon . bx,by - type of ghost nodes the array have. 809bff4a2f0SMatthew G. Knepley Use one of DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_PERIODIC. 810aa219208SBarry Smith . stencil_type - stencil type. Use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR. 81147c6ae99SBarry Smith . M,N - global dimension in each direction of the array (use -M and or -N to indicate that it may be set to a different value 81247c6ae99SBarry Smith from the command line with -da_grid_x <M> -da_grid_y <N>) 81347c6ae99SBarry Smith . m,n - corresponding number of processors in each dimension 81447c6ae99SBarry Smith (or PETSC_DECIDE to have calculated) 81547c6ae99SBarry Smith . dof - number of degrees of freedom per node 81647c6ae99SBarry Smith . s - stencil width 81747c6ae99SBarry Smith - lx, ly - arrays containing the number of nodes in each cell along 8180298fd71SBarry Smith the x and y coordinates, or NULL. If non-null, these 81947c6ae99SBarry Smith must be of length as m and n, and the corresponding 82047c6ae99SBarry Smith m and n cannot be PETSC_DECIDE. The sum of the lx[] entries 82147c6ae99SBarry Smith must be M, and the sum of the ly[] entries must be N. 82247c6ae99SBarry Smith 82347c6ae99SBarry Smith Output Parameter: 82447c6ae99SBarry Smith . da - the resulting distributed array object 82547c6ae99SBarry Smith 82647c6ae99SBarry Smith Options Database Key: 827706a11cbSBarry Smith + -dm_view - Calls DMView() at the conclusion of DMDACreate2d() 82847c6ae99SBarry Smith . -da_grid_x <nx> - number of grid points in x direction, if M < 0 82947c6ae99SBarry Smith . -da_grid_y <ny> - number of grid points in y direction, if N < 0 83047c6ae99SBarry Smith . -da_processors_x <nx> - number of processors in x direction 83147c6ae99SBarry Smith . -da_processors_y <ny> - number of processors in y direction 832e0f5d30fSBarry Smith . -da_refine_x <rx> - refinement ratio in x direction 833e0f5d30fSBarry Smith . -da_refine_y <ry> - refinement ratio in y direction 834e0f5d30fSBarry Smith - -da_refine <n> - refine the DMDA n times before creating, if M or N < 0 835e0f5d30fSBarry Smith 83647c6ae99SBarry Smith 83747c6ae99SBarry Smith Level: beginner 83847c6ae99SBarry Smith 83947c6ae99SBarry Smith Notes: 840aa219208SBarry Smith The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the 841aa219208SBarry Smith standard 5-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes 84247c6ae99SBarry Smith the standard 9-pt stencil. 84347c6ae99SBarry Smith 844aa219208SBarry Smith The array data itself is NOT stored in the DMDA, it is stored in Vec objects; 845564755cdSBarry Smith The appropriate vector objects can be obtained with calls to DMCreateGlobalVector() 846564755cdSBarry Smith and DMCreateLocalVector() and calls to VecDuplicate() if more are needed. 84747c6ae99SBarry 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); 863aa219208SBarry Smith ierr = DMDASetDim(*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 /* This violates the behavior for other classes, but right now users expect negative dimensions to be handled this way */ 8729a42bb27SBarry Smith ierr = DMSetFromOptions(*da);CHKERRQ(ierr); 8739a42bb27SBarry Smith ierr = DMSetUp(*da);CHKERRQ(ierr); 87447c6ae99SBarry Smith PetscFunctionReturn(0); 87547c6ae99SBarry Smith } 876