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; 9*c9493c35SStefano 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); 20*c9493c35SStefano 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); 298135c375SStefano Zampini if (format != PETSC_VIEWER_ASCII_VTK && format != PETSC_VIEWER_ASCII_VTK_CELL && format != PETSC_VIEWER_ASCII_GLVIS) { 30aa219208SBarry Smith DMDALocalInfo info; 31aa219208SBarry Smith ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr); 321575c14dSBarry Smith ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 3347c6ae99SBarry 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); 3447c6ae99SBarry 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); 3547c6ae99SBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 361575c14dSBarry Smith ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 378135c375SStefano Zampini } else if (format == PETSC_VIEWER_ASCII_GLVIS) { 388135c375SStefano Zampini ierr = DMView_DA_GLVis(da,viewer);CHKERRQ(ierr); 393da9ae13SJed Brown } else { 403da9ae13SJed Brown ierr = DMView_DA_VTK(da,viewer);CHKERRQ(ierr); 4147c6ae99SBarry Smith } 4247c6ae99SBarry Smith } else if (isdraw) { 4347c6ae99SBarry Smith PetscDraw draw; 4447c6ae99SBarry Smith double ymin = -1*dd->s-1,ymax = dd->N+dd->s; 4547c6ae99SBarry Smith double xmin = -1*dd->s-1,xmax = dd->M+dd->s; 4647c6ae99SBarry Smith double x,y; 478ea3bf28SBarry Smith PetscInt base; 488ea3bf28SBarry Smith const PetscInt *idx; 4947c6ae99SBarry Smith char node[10]; 5047c6ae99SBarry Smith PetscBool isnull; 5147c6ae99SBarry Smith 5247c6ae99SBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 535b399a63SLisandro Dalcin ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 545b399a63SLisandro Dalcin if (isnull) PetscFunctionReturn(0); 5547c6ae99SBarry Smith 565b399a63SLisandro Dalcin ierr = PetscDrawCheckResizedWindow(draw);CHKERRQ(ierr); 575b399a63SLisandro Dalcin ierr = PetscDrawClear(draw);CHKERRQ(ierr); 585b399a63SLisandro Dalcin ierr = PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);CHKERRQ(ierr); 595b399a63SLisandro Dalcin 605b399a63SLisandro Dalcin ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr); 6147c6ae99SBarry Smith /* first processor draw all node lines */ 6247c6ae99SBarry Smith if (!rank) { 6347c6ae99SBarry Smith ymin = 0.0; ymax = dd->N - 1; 6447c6ae99SBarry Smith for (xmin=0; xmin<dd->M; xmin++) { 6547c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr); 6647c6ae99SBarry Smith } 6747c6ae99SBarry Smith xmin = 0.0; xmax = dd->M - 1; 6847c6ae99SBarry Smith for (ymin=0; ymin<dd->N; ymin++) { 6947c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);CHKERRQ(ierr); 7047c6ae99SBarry Smith } 7147c6ae99SBarry Smith } 725b399a63SLisandro Dalcin ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr); 735b399a63SLisandro Dalcin ierr = PetscDrawFlush(draw);CHKERRQ(ierr); 7447c6ae99SBarry Smith ierr = PetscDrawPause(draw);CHKERRQ(ierr); 7547c6ae99SBarry Smith 765b399a63SLisandro Dalcin ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr); 7747c6ae99SBarry Smith /* draw my box */ 785b399a63SLisandro Dalcin xmin = dd->xs/dd->w; xmax =(dd->xe-1)/dd->w; ymin = dd->ys; ymax = dd->ye - 1; 7947c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);CHKERRQ(ierr); 8047c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 8147c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 8247c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 8347c6ae99SBarry Smith /* put in numbers */ 8447c6ae99SBarry Smith base = (dd->base)/dd->w; 8547c6ae99SBarry Smith for (y=ymin; y<=ymax; y++) { 8647c6ae99SBarry Smith for (x=xmin; x<=xmax; x++) { 875b399a63SLisandro Dalcin ierr = PetscSNPrintf(node,sizeof(node),"%d",(int)base++);CHKERRQ(ierr); 8847c6ae99SBarry Smith ierr = PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node);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); 965b399a63SLisandro Dalcin /* overlay ghost numbers, useful for error checking */ 9745b6f7e9SBarry Smith ierr = ISLocalToGlobalMappingGetBlockIndices(da->ltogmap,&idx);CHKERRQ(ierr); 985b399a63SLisandro Dalcin base = 0; xmin = dd->Xs; xmax = dd->Xe; ymin = dd->Ys; ymax = dd->Ye; 9947c6ae99SBarry Smith for (y=ymin; y<ymax; y++) { 10047c6ae99SBarry Smith for (x=xmin; x<xmax; x++) { 10147c6ae99SBarry Smith if ((base % dd->w) == 0) { 1025b399a63SLisandro Dalcin ierr = PetscSNPrintf(node,sizeof(node),"%d",(int)(idx[base/dd->w]));CHKERRQ(ierr); 10347c6ae99SBarry Smith ierr = PetscDrawString(draw,x/dd->w,y,PETSC_DRAW_BLUE,node);CHKERRQ(ierr); 10447c6ae99SBarry Smith } 10547c6ae99SBarry Smith base++; 10647c6ae99SBarry Smith } 10747c6ae99SBarry Smith } 108302440fdSBarry Smith ierr = ISLocalToGlobalMappingRestoreBlockIndices(da->ltogmap,&idx);CHKERRQ(ierr); 1095b399a63SLisandro Dalcin ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr); 1105b399a63SLisandro Dalcin ierr = PetscDrawFlush(draw);CHKERRQ(ierr); 11147c6ae99SBarry Smith ierr = PetscDrawPause(draw);CHKERRQ(ierr); 112832b7cebSLisandro Dalcin ierr = PetscDrawSave(draw);CHKERRQ(ierr); 113*c9493c35SStefano Zampini } else if (isglvis) { 114*c9493c35SStefano Zampini ierr = DMView_DA_GLVis(da,viewer);CHKERRQ(ierr); 1159a42bb27SBarry Smith } else if (isbinary) { 1169a42bb27SBarry Smith ierr = DMView_DA_Binary(da,viewer);CHKERRQ(ierr); 1179a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 1189a42bb27SBarry Smith } else if (ismatlab) { 1199a42bb27SBarry Smith ierr = DMView_DA_Matlab(da,viewer);CHKERRQ(ierr); 1209a42bb27SBarry Smith #endif 12111aeaf0aSBarry Smith } 12247c6ae99SBarry Smith PetscFunctionReturn(0); 12347c6ae99SBarry Smith } 12447c6ae99SBarry Smith 12547c6ae99SBarry Smith /* 12647c6ae99SBarry Smith M is number of grid points 12747c6ae99SBarry Smith m is number of processors 12847c6ae99SBarry Smith 12947c6ae99SBarry Smith */ 1307087cfbeSBarry Smith PetscErrorCode DMDASplitComm2d(MPI_Comm comm,PetscInt M,PetscInt N,PetscInt sw,MPI_Comm *outcomm) 13147c6ae99SBarry Smith { 13247c6ae99SBarry Smith PetscErrorCode ierr; 13347c6ae99SBarry Smith PetscInt m,n = 0,x = 0,y = 0; 13447c6ae99SBarry Smith PetscMPIInt size,csize,rank; 13547c6ae99SBarry Smith 13647c6ae99SBarry Smith PetscFunctionBegin; 13747c6ae99SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 13847c6ae99SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 13947c6ae99SBarry Smith 14047c6ae99SBarry Smith csize = 4*size; 14147c6ae99SBarry Smith do { 14247c6ae99SBarry 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); 14347c6ae99SBarry Smith csize = csize/4; 14447c6ae99SBarry Smith 145369cc0aeSBarry Smith m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)csize)/((PetscReal)N))); 14647c6ae99SBarry Smith if (!m) m = 1; 14747c6ae99SBarry Smith while (m > 0) { 14847c6ae99SBarry Smith n = csize/m; 14947c6ae99SBarry Smith if (m*n == csize) break; 15047c6ae99SBarry Smith m--; 15147c6ae99SBarry Smith } 15247c6ae99SBarry Smith if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;} 15347c6ae99SBarry Smith 15447c6ae99SBarry Smith x = M/m + ((M % m) > ((csize-1) % m)); 15547c6ae99SBarry Smith y = (N + (csize-1)/m)/n; 15647c6ae99SBarry Smith } while ((x < 4 || y < 4) && csize > 1); 15747c6ae99SBarry Smith if (size != csize) { 15847c6ae99SBarry Smith MPI_Group entire_group,sub_group; 15947c6ae99SBarry Smith PetscMPIInt i,*groupies; 16047c6ae99SBarry Smith 16147c6ae99SBarry Smith ierr = MPI_Comm_group(comm,&entire_group);CHKERRQ(ierr); 162785e854fSJed Brown ierr = PetscMalloc1(csize,&groupies);CHKERRQ(ierr); 16347c6ae99SBarry Smith for (i=0; i<csize; i++) { 16447c6ae99SBarry Smith groupies[i] = (rank/csize)*csize + i; 16547c6ae99SBarry Smith } 16647c6ae99SBarry Smith ierr = MPI_Group_incl(entire_group,csize,groupies,&sub_group);CHKERRQ(ierr); 16747c6ae99SBarry Smith ierr = PetscFree(groupies);CHKERRQ(ierr); 16847c6ae99SBarry Smith ierr = MPI_Comm_create(comm,sub_group,outcomm);CHKERRQ(ierr); 16947c6ae99SBarry Smith ierr = MPI_Group_free(&entire_group);CHKERRQ(ierr); 17047c6ae99SBarry Smith ierr = MPI_Group_free(&sub_group);CHKERRQ(ierr); 171aa219208SBarry Smith ierr = PetscInfo1(0,"DMDASplitComm2d:Creating redundant coarse problems of size %d\n",csize);CHKERRQ(ierr); 17247c6ae99SBarry Smith } else { 17347c6ae99SBarry Smith *outcomm = comm; 17447c6ae99SBarry Smith } 17547c6ae99SBarry Smith PetscFunctionReturn(0); 17647c6ae99SBarry Smith } 17747c6ae99SBarry Smith 17847c6ae99SBarry Smith #if defined(new) 17947c6ae99SBarry Smith /* 180aa219208SBarry Smith DMDAGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix where local 181aa219208SBarry Smith function lives on a DMDA 18247c6ae99SBarry Smith 18347c6ae99SBarry Smith y ~= (F(u + ha) - F(u))/h, 18447c6ae99SBarry Smith where F = nonlinear function, as set by SNESSetFunction() 18547c6ae99SBarry Smith u = current iterate 18647c6ae99SBarry Smith h = difference interval 18747c6ae99SBarry Smith */ 188aa219208SBarry Smith PetscErrorCode DMDAGetDiagonal_MFFD(DM da,Vec U,Vec a) 18947c6ae99SBarry Smith { 19047c6ae99SBarry Smith PetscScalar h,*aa,*ww,v; 19147c6ae99SBarry Smith PetscReal epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON; 19247c6ae99SBarry Smith PetscErrorCode ierr; 19347c6ae99SBarry Smith PetscInt gI,nI; 19447c6ae99SBarry Smith MatStencil stencil; 195aa219208SBarry Smith DMDALocalInfo info; 19647c6ae99SBarry Smith 19747c6ae99SBarry Smith PetscFunctionBegin; 19847c6ae99SBarry Smith ierr = (*ctx->func)(0,U,a,ctx->funcctx);CHKERRQ(ierr); 19947c6ae99SBarry Smith ierr = (*ctx->funcisetbase)(U,ctx->funcctx);CHKERRQ(ierr); 20047c6ae99SBarry Smith 20147c6ae99SBarry Smith ierr = VecGetArray(U,&ww);CHKERRQ(ierr); 20247c6ae99SBarry Smith ierr = VecGetArray(a,&aa);CHKERRQ(ierr); 20347c6ae99SBarry Smith 20447c6ae99SBarry Smith nI = 0; 20547c6ae99SBarry Smith h = ww[gI]; 20647c6ae99SBarry Smith if (h == 0.0) h = 1.0; 20747c6ae99SBarry Smith if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin; 20847c6ae99SBarry Smith else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin; 20947c6ae99SBarry Smith h *= epsilon; 21047c6ae99SBarry Smith 21147c6ae99SBarry Smith ww[gI] += h; 21247c6ae99SBarry Smith ierr = (*ctx->funci)(i,w,&v,ctx->funcctx);CHKERRQ(ierr); 21347c6ae99SBarry Smith aa[nI] = (v - aa[nI])/h; 21447c6ae99SBarry Smith ww[gI] -= h; 21547c6ae99SBarry Smith nI++; 2168865f1eaSKarl Rupp 21747c6ae99SBarry Smith ierr = VecRestoreArray(U,&ww);CHKERRQ(ierr); 21847c6ae99SBarry Smith ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr); 21947c6ae99SBarry Smith PetscFunctionReturn(0); 22047c6ae99SBarry Smith } 22147c6ae99SBarry Smith #endif 22247c6ae99SBarry Smith 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; 239bd1fc5aeSBarry Smith PetscInt xs,xe,ys,ye,x,y,Xs,Xe,Ys,Ye,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; 245bd1fc5aeSBarry Smith VecScatter 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) 253bafee8b4SSatish 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); 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 2627d310018SBarry Smith dd->p = 1; 26347c6ae99SBarry Smith if (m != PETSC_DECIDE) { 26447c6ae99SBarry Smith if (m < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m); 26547c6ae99SBarry Smith else if (m > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size); 26647c6ae99SBarry Smith } 26747c6ae99SBarry Smith if (n != PETSC_DECIDE) { 26847c6ae99SBarry Smith if (n < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n); 26947c6ae99SBarry Smith else if (n > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size); 27047c6ae99SBarry Smith } 27147c6ae99SBarry Smith 27247c6ae99SBarry Smith if (m == PETSC_DECIDE || n == PETSC_DECIDE) { 27347c6ae99SBarry Smith if (n != PETSC_DECIDE) { 27447c6ae99SBarry Smith m = size/n; 27547c6ae99SBarry Smith } else if (m != PETSC_DECIDE) { 27647c6ae99SBarry Smith n = size/m; 27747c6ae99SBarry Smith } else { 27847c6ae99SBarry Smith /* try for squarish distribution */ 279369cc0aeSBarry Smith m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)N))); 28047c6ae99SBarry Smith if (!m) m = 1; 28147c6ae99SBarry Smith while (m > 0) { 28247c6ae99SBarry Smith n = size/m; 28347c6ae99SBarry Smith if (m*n == size) break; 28447c6ae99SBarry Smith m--; 28547c6ae99SBarry Smith } 28647c6ae99SBarry Smith if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;} 28747c6ae99SBarry Smith } 28847c6ae99SBarry 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 "); 28947c6ae99SBarry Smith } else if (m*n != size) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition"); 29047c6ae99SBarry Smith 29147c6ae99SBarry Smith if (M < m) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m); 29247c6ae99SBarry Smith if (N < n) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n); 29347c6ae99SBarry Smith 29447c6ae99SBarry Smith /* 29547c6ae99SBarry Smith Determine locally owned region 29647c6ae99SBarry Smith xs is the first local node number, x is the number of local nodes 29747c6ae99SBarry Smith */ 29847c6ae99SBarry Smith if (!lx) { 299785e854fSJed Brown ierr = PetscMalloc1(m, &dd->lx);CHKERRQ(ierr); 30047c6ae99SBarry Smith lx = dd->lx; 30147c6ae99SBarry Smith for (i=0; i<m; i++) { 30247c6ae99SBarry Smith lx[i] = M/m + ((M % m) > i); 30347c6ae99SBarry Smith } 30447c6ae99SBarry Smith } 30547c6ae99SBarry Smith x = lx[rank % m]; 30647c6ae99SBarry Smith xs = 0; 30747c6ae99SBarry Smith for (i=0; i<(rank % m); i++) { 30847c6ae99SBarry Smith xs += lx[i]; 30947c6ae99SBarry Smith } 31047c6ae99SBarry Smith #if defined(PETSC_USE_DEBUG) 31147c6ae99SBarry Smith left = xs; 31247c6ae99SBarry Smith for (i=(rank % m); i<m; i++) { 31347c6ae99SBarry Smith left += lx[i]; 31447c6ae99SBarry Smith } 31547c6ae99SBarry 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); 31647c6ae99SBarry Smith #endif 31747c6ae99SBarry Smith 31847c6ae99SBarry Smith /* 31947c6ae99SBarry Smith Determine locally owned region 32047c6ae99SBarry Smith ys is the first local node number, y is the number of local nodes 32147c6ae99SBarry Smith */ 32247c6ae99SBarry Smith if (!ly) { 323785e854fSJed Brown ierr = PetscMalloc1(n, &dd->ly);CHKERRQ(ierr); 32447c6ae99SBarry Smith ly = dd->ly; 32547c6ae99SBarry Smith for (i=0; i<n; i++) { 32647c6ae99SBarry Smith ly[i] = N/n + ((N % n) > i); 32747c6ae99SBarry Smith } 32847c6ae99SBarry Smith } 32947c6ae99SBarry Smith y = ly[rank/m]; 33047c6ae99SBarry Smith ys = 0; 33147c6ae99SBarry Smith for (i=0; i<(rank/m); i++) { 33247c6ae99SBarry Smith ys += ly[i]; 33347c6ae99SBarry Smith } 33447c6ae99SBarry Smith #if defined(PETSC_USE_DEBUG) 33547c6ae99SBarry Smith left = ys; 33647c6ae99SBarry Smith for (i=(rank/m); i<n; i++) { 33747c6ae99SBarry Smith left += ly[i]; 33847c6ae99SBarry Smith } 33947c6ae99SBarry 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); 34047c6ae99SBarry Smith #endif 34147c6ae99SBarry Smith 342bcea557cSEthan Coon /* 343bcea557cSEthan Coon check if the scatter requires more than one process neighbor or wraps around 344bcea557cSEthan Coon the domain more than once 345bcea557cSEthan Coon */ 346bff4a2f0SMatthew 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); 347bff4a2f0SMatthew 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); 34847c6ae99SBarry Smith xe = xs + x; 34947c6ae99SBarry Smith ye = ys + y; 35047c6ae99SBarry Smith 351ce00eea3SSatish Balay /* determine ghost region (Xs) and region scattered into (IXs) */ 352d9c9ebe5SPeter Brune if (xs-s > 0) { 353d9c9ebe5SPeter Brune Xs = xs - s; IXs = xs - s; 35488661749SPeter Brune } else { 35588661749SPeter Brune if (bx) { 35688661749SPeter Brune Xs = xs - s; 35788661749SPeter Brune } else { 35888661749SPeter Brune Xs = 0; 35988661749SPeter Brune } 36088661749SPeter Brune IXs = 0; 36188661749SPeter Brune } 362d9c9ebe5SPeter Brune if (xe+s <= M) { 363d9c9ebe5SPeter Brune Xe = xe + s; IXe = xe + s; 36488661749SPeter Brune } else { 36588661749SPeter Brune if (bx) { 366d9c9ebe5SPeter Brune Xs = xs - s; Xe = xe + s; 36788661749SPeter Brune } else { 36888661749SPeter Brune Xe = M; 36988661749SPeter Brune } 37088661749SPeter Brune IXe = M; 37188661749SPeter Brune } 37247c6ae99SBarry Smith 373bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC || bx == DM_BOUNDARY_MIRROR) { 374d9c9ebe5SPeter Brune IXs = xs - s; 375d9c9ebe5SPeter Brune IXe = xe + s; 376d9c9ebe5SPeter Brune Xs = xs - s; 377d9c9ebe5SPeter Brune Xe = xe + s; 37888661749SPeter Brune } 37947c6ae99SBarry Smith 380d9c9ebe5SPeter Brune if (ys-s > 0) { 381d9c9ebe5SPeter Brune Ys = ys - s; IYs = ys - s; 38288661749SPeter Brune } else { 38388661749SPeter Brune if (by) { 38488661749SPeter Brune Ys = ys - s; 38588661749SPeter Brune } else { 38688661749SPeter Brune Ys = 0; 38788661749SPeter Brune } 38888661749SPeter Brune IYs = 0; 38988661749SPeter Brune } 390d9c9ebe5SPeter Brune if (ye+s <= N) { 391d9c9ebe5SPeter Brune Ye = ye + s; IYe = ye + s; 39288661749SPeter Brune } else { 39388661749SPeter Brune if (by) { 39488661749SPeter Brune Ye = ye + s; 39588661749SPeter Brune } else { 39688661749SPeter Brune Ye = N; 39788661749SPeter Brune } 39888661749SPeter Brune IYe = N; 39988661749SPeter Brune } 40088661749SPeter Brune 401bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_PERIODIC || by == DM_BOUNDARY_MIRROR) { 402d9c9ebe5SPeter Brune IYs = ys - s; 403d9c9ebe5SPeter Brune IYe = ye + s; 404d9c9ebe5SPeter Brune Ys = ys - s; 405d9c9ebe5SPeter Brune Ye = ye + s; 40688661749SPeter Brune } 40788661749SPeter Brune 40888661749SPeter Brune /* stencil length in each direction */ 409d9c9ebe5SPeter Brune s_x = s; 410d9c9ebe5SPeter Brune s_y = s; 41147c6ae99SBarry Smith 41247c6ae99SBarry Smith /* determine starting point of each processor */ 41347c6ae99SBarry Smith nn = x*y; 414dcca6d9dSJed Brown ierr = PetscMalloc2(size+1,&bases,size,&ldims);CHKERRQ(ierr); 41547c6ae99SBarry Smith ierr = MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);CHKERRQ(ierr); 41647c6ae99SBarry Smith bases[0] = 0; 41747c6ae99SBarry Smith for (i=1; i<=size; i++) { 41847c6ae99SBarry Smith bases[i] = ldims[i-1]; 41947c6ae99SBarry Smith } 42047c6ae99SBarry Smith for (i=1; i<=size; i++) { 42147c6ae99SBarry Smith bases[i] += bases[i-1]; 42247c6ae99SBarry Smith } 423ce00eea3SSatish Balay base = bases[rank]*dof; 42447c6ae99SBarry Smith 42547c6ae99SBarry Smith /* allocate the base parallel and sequential vectors */ 426ce00eea3SSatish Balay dd->Nlocal = x*y*dof; 427b1fb7eb7SBarry Smith ierr = VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,NULL,&global);CHKERRQ(ierr); 428ce00eea3SSatish Balay dd->nlocal = (Xe-Xs)*(Ye-Ys)*dof; 429b1fb7eb7SBarry Smith ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,NULL,&local);CHKERRQ(ierr); 43047c6ae99SBarry Smith 43147c6ae99SBarry Smith /* generate appropriate vector scatters */ 43247c6ae99SBarry Smith /* local to global inserts non-ghost point region into global */ 43302fe608eSBarry Smith ierr = PetscMalloc1((IXe-IXs)*(IYe-IYs),&idx);CHKERRQ(ierr); 434ce00eea3SSatish Balay left = xs - Xs; right = left + x; 435ce00eea3SSatish Balay down = ys - Ys; up = down + y; 43647c6ae99SBarry Smith count = 0; 43747c6ae99SBarry Smith for (i=down; i<up; i++) { 438ce00eea3SSatish Balay for (j=left; j<right; j++) { 439ce00eea3SSatish Balay idx[count++] = i*(Xe-Xs) + j; 44047c6ae99SBarry Smith } 44147c6ae99SBarry Smith } 44247c6ae99SBarry Smith 443ce00eea3SSatish Balay /* global to local must include ghost points within the domain, 444ce00eea3SSatish Balay but not ghost points outside the domain that aren't periodic */ 445d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_BOX) { 446ce00eea3SSatish Balay left = IXs - Xs; right = left + (IXe-IXs); 447ce00eea3SSatish Balay down = IYs - Ys; up = down + (IYe-IYs); 448ce00eea3SSatish Balay count = 0; 449ce00eea3SSatish Balay for (i=down; i<up; i++) { 450ce00eea3SSatish Balay for (j=left; j<right; j++) { 451ce00eea3SSatish Balay idx[count++] = j + i*(Xe-Xs); 452ce00eea3SSatish Balay } 453ce00eea3SSatish Balay } 454ce00eea3SSatish Balay ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr); 455ce00eea3SSatish Balay 45647c6ae99SBarry Smith } else { 45747c6ae99SBarry Smith /* must drop into cross shape region */ 45847c6ae99SBarry Smith /* ---------| 45947c6ae99SBarry Smith | top | 460ce00eea3SSatish Balay |--- ---| up 46147c6ae99SBarry Smith | middle | 46247c6ae99SBarry Smith | | 463ce00eea3SSatish Balay ---- ---- down 46447c6ae99SBarry Smith | bottom | 46547c6ae99SBarry Smith ----------- 46647c6ae99SBarry Smith Xs xs xe Xe */ 467ce00eea3SSatish Balay left = xs - Xs; right = left + x; 468ce00eea3SSatish Balay down = ys - Ys; up = down + y; 46947c6ae99SBarry Smith count = 0; 470ce00eea3SSatish Balay /* bottom */ 471ce00eea3SSatish Balay for (i=(IYs-Ys); i<down; i++) { 472ce00eea3SSatish Balay for (j=left; j<right; j++) { 473ce00eea3SSatish Balay idx[count++] = j + i*(Xe-Xs); 47447c6ae99SBarry Smith } 47547c6ae99SBarry Smith } 47647c6ae99SBarry Smith /* middle */ 47747c6ae99SBarry Smith for (i=down; i<up; i++) { 478ce00eea3SSatish Balay for (j=(IXs-Xs); j<(IXe-Xs); j++) { 479ce00eea3SSatish Balay idx[count++] = j + i*(Xe-Xs); 48047c6ae99SBarry Smith } 48147c6ae99SBarry Smith } 48247c6ae99SBarry Smith /* top */ 483ce00eea3SSatish Balay for (i=up; i<up+IYe-ye; i++) { 484ce00eea3SSatish Balay for (j=left; j<right; j++) { 485ce00eea3SSatish Balay idx[count++] = j + i*(Xe-Xs); 48647c6ae99SBarry Smith } 48747c6ae99SBarry Smith } 48847c6ae99SBarry Smith ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr); 48947c6ae99SBarry Smith } 49047c6ae99SBarry Smith 49147c6ae99SBarry Smith 49247c6ae99SBarry Smith /* determine who lies on each side of us stored in n6 n7 n8 49347c6ae99SBarry Smith n3 n5 49447c6ae99SBarry Smith n0 n1 n2 49547c6ae99SBarry Smith */ 49647c6ae99SBarry Smith 49747c6ae99SBarry Smith /* Assume the Non-Periodic Case */ 49847c6ae99SBarry Smith n1 = rank - m; 49947c6ae99SBarry Smith if (rank % m) { 50047c6ae99SBarry Smith n0 = n1 - 1; 50147c6ae99SBarry Smith } else { 50247c6ae99SBarry Smith n0 = -1; 50347c6ae99SBarry Smith } 50447c6ae99SBarry Smith if ((rank+1) % m) { 50547c6ae99SBarry Smith n2 = n1 + 1; 50647c6ae99SBarry Smith n5 = rank + 1; 50747c6ae99SBarry Smith n8 = rank + m + 1; if (n8 >= m*n) n8 = -1; 50847c6ae99SBarry Smith } else { 50947c6ae99SBarry Smith n2 = -1; n5 = -1; n8 = -1; 51047c6ae99SBarry Smith } 51147c6ae99SBarry Smith if (rank % m) { 51247c6ae99SBarry Smith n3 = rank - 1; 51347c6ae99SBarry Smith n6 = n3 + m; if (n6 >= m*n) n6 = -1; 51447c6ae99SBarry Smith } else { 51547c6ae99SBarry Smith n3 = -1; n6 = -1; 51647c6ae99SBarry Smith } 51747c6ae99SBarry Smith n7 = rank + m; if (n7 >= m*n) n7 = -1; 51847c6ae99SBarry Smith 519bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC && by == DM_BOUNDARY_PERIODIC) { 52047c6ae99SBarry Smith /* Modify for Periodic Cases */ 52147c6ae99SBarry Smith /* Handle all four corners */ 52247c6ae99SBarry Smith if ((n6 < 0) && (n7 < 0) && (n3 < 0)) n6 = m-1; 52347c6ae99SBarry Smith if ((n8 < 0) && (n7 < 0) && (n5 < 0)) n8 = 0; 52447c6ae99SBarry Smith if ((n2 < 0) && (n5 < 0) && (n1 < 0)) n2 = size-m; 52547c6ae99SBarry Smith if ((n0 < 0) && (n3 < 0) && (n1 < 0)) n0 = size-1; 52647c6ae99SBarry Smith 52747c6ae99SBarry Smith /* Handle Top and Bottom Sides */ 52847c6ae99SBarry Smith if (n1 < 0) n1 = rank + m * (n-1); 52947c6ae99SBarry Smith if (n7 < 0) n7 = rank - m * (n-1); 53047c6ae99SBarry Smith if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1; 53147c6ae99SBarry Smith if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1; 53247c6ae99SBarry Smith if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1; 53347c6ae99SBarry Smith if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1; 53447c6ae99SBarry Smith 53547c6ae99SBarry Smith /* Handle Left and Right Sides */ 53647c6ae99SBarry Smith if (n3 < 0) n3 = rank + (m-1); 53747c6ae99SBarry Smith if (n5 < 0) n5 = rank - (m-1); 53847c6ae99SBarry Smith if ((n1 >= 0) && (n0 < 0)) n0 = rank-1; 53947c6ae99SBarry Smith if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1; 54047c6ae99SBarry Smith if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1; 54147c6ae99SBarry Smith if ((n7 >= 0) && (n8 < 0)) n8 = rank+1; 542bff4a2f0SMatthew G. Knepley } else if (by == DM_BOUNDARY_PERIODIC) { /* Handle Top and Bottom Sides */ 543ce00eea3SSatish Balay if (n1 < 0) n1 = rank + m * (n-1); 544ce00eea3SSatish Balay if (n7 < 0) n7 = rank - m * (n-1); 545ce00eea3SSatish Balay if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1; 546ce00eea3SSatish Balay if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1; 547ce00eea3SSatish Balay if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1; 548ce00eea3SSatish Balay if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1; 549bff4a2f0SMatthew G. Knepley } else if (bx == DM_BOUNDARY_PERIODIC) { /* Handle Left and Right Sides */ 550ce00eea3SSatish Balay if (n3 < 0) n3 = rank + (m-1); 551ce00eea3SSatish Balay if (n5 < 0) n5 = rank - (m-1); 552ce00eea3SSatish Balay if ((n1 >= 0) && (n0 < 0)) n0 = rank-1; 553ce00eea3SSatish Balay if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1; 554ce00eea3SSatish Balay if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1; 555ce00eea3SSatish Balay if ((n7 >= 0) && (n8 < 0)) n8 = rank+1; 55647c6ae99SBarry Smith } 557ce00eea3SSatish Balay 558785e854fSJed Brown ierr = PetscMalloc1(9,&dd->neighbors);CHKERRQ(ierr); 5598865f1eaSKarl Rupp 56047c6ae99SBarry Smith dd->neighbors[0] = n0; 56147c6ae99SBarry Smith dd->neighbors[1] = n1; 56247c6ae99SBarry Smith dd->neighbors[2] = n2; 56347c6ae99SBarry Smith dd->neighbors[3] = n3; 56447c6ae99SBarry Smith dd->neighbors[4] = rank; 56547c6ae99SBarry Smith dd->neighbors[5] = n5; 56647c6ae99SBarry Smith dd->neighbors[6] = n6; 56747c6ae99SBarry Smith dd->neighbors[7] = n7; 56847c6ae99SBarry Smith dd->neighbors[8] = n8; 56947c6ae99SBarry Smith 570d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_STAR) { 57147c6ae99SBarry Smith /* save corner processor numbers */ 57247c6ae99SBarry Smith sn0 = n0; sn2 = n2; sn6 = n6; sn8 = n8; 57347c6ae99SBarry Smith n0 = n2 = n6 = n8 = -1; 57447c6ae99SBarry Smith } 57547c6ae99SBarry Smith 576785e854fSJed Brown ierr = PetscMalloc1((Xe-Xs)*(Ye-Ys),&idx);CHKERRQ(ierr); 57747c6ae99SBarry Smith 578ce00eea3SSatish Balay nn = 0; 57947c6ae99SBarry Smith xbase = bases[rank]; 58047c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 58147c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 582ce00eea3SSatish Balay x_t = lx[n0 % m]; 58347c6ae99SBarry Smith y_t = ly[(n0/m)]; 58447c6ae99SBarry Smith s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x; 5858865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 58647c6ae99SBarry Smith } 587ac119b13SBarry Smith 58847c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 58947c6ae99SBarry Smith x_t = x; 59047c6ae99SBarry Smith y_t = ly[(n1/m)]; 59147c6ae99SBarry Smith s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t; 5928865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 593bff4a2f0SMatthew G. Knepley } else if (by == DM_BOUNDARY_MIRROR) { 5948865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(s_y - i + 1) + j; 59547c6ae99SBarry Smith } 596ac119b13SBarry Smith 59747c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 598ce00eea3SSatish Balay x_t = lx[n2 % m]; 59947c6ae99SBarry Smith y_t = ly[(n2/m)]; 60047c6ae99SBarry Smith s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t; 6018865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 60247c6ae99SBarry Smith } 60347c6ae99SBarry Smith } 60447c6ae99SBarry Smith 60547c6ae99SBarry Smith for (i=0; i<y; i++) { 60647c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 607ce00eea3SSatish Balay x_t = lx[n3 % m]; 60847c6ae99SBarry Smith /* y_t = y; */ 60947c6ae99SBarry Smith s_t = bases[n3] + (i+1)*x_t - s_x; 6108865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 611bff4a2f0SMatthew G. Knepley } else if (bx == DM_BOUNDARY_MIRROR) { 6128865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*i + s_x - j; 61347c6ae99SBarry Smith } 61447c6ae99SBarry Smith 6158865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = xbase++; /* interior */ 61647c6ae99SBarry Smith 61747c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 618ce00eea3SSatish Balay x_t = lx[n5 % m]; 61947c6ae99SBarry Smith /* y_t = y; */ 62047c6ae99SBarry Smith s_t = bases[n5] + (i)*x_t; 6218865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 622bff4a2f0SMatthew G. Knepley } else if (bx == DM_BOUNDARY_MIRROR) { 6238865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*(i + 1) - 2 - j; 62447c6ae99SBarry Smith } 62547c6ae99SBarry Smith } 62647c6ae99SBarry Smith 62747c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 62847c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 629ce00eea3SSatish Balay x_t = lx[n6 % m]; 63047c6ae99SBarry Smith /* y_t = ly[(n6/m)]; */ 63147c6ae99SBarry Smith s_t = bases[n6] + (i)*x_t - s_x; 6328865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 63347c6ae99SBarry Smith } 634ac119b13SBarry Smith 63547c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 63647c6ae99SBarry Smith x_t = x; 63747c6ae99SBarry Smith /* y_t = ly[(n7/m)]; */ 63847c6ae99SBarry Smith s_t = bases[n7] + (i-1)*x_t; 6398865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 640bff4a2f0SMatthew G. Knepley } else if (by == DM_BOUNDARY_MIRROR) { 6418865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(y - i - 1) + j; 64247c6ae99SBarry Smith } 643ac119b13SBarry Smith 64447c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 645ce00eea3SSatish Balay x_t = lx[n8 % m]; 64647c6ae99SBarry Smith /* y_t = ly[(n8/m)]; */ 64747c6ae99SBarry Smith s_t = bases[n8] + (i-1)*x_t; 6488865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 64947c6ae99SBarry Smith } 65047c6ae99SBarry Smith } 65147c6ae99SBarry Smith 652b1fb7eb7SBarry Smith ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_USE_POINTER,&from);CHKERRQ(ierr); 65347c6ae99SBarry Smith ierr = VecScatterCreate(global,from,local,to,>ol);CHKERRQ(ierr); 6543bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)gtol);CHKERRQ(ierr); 655fcfd50ebSBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 656fcfd50ebSBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 65747c6ae99SBarry Smith 658d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_STAR) { 659ce00eea3SSatish Balay n0 = sn0; n2 = sn2; n6 = sn6; n8 = sn8; 660ce00eea3SSatish Balay } 661ce00eea3SSatish Balay 662288e7d53SBarry Smith if (((stencil_type == DMDA_STENCIL_STAR) || (bx && bx != DM_BOUNDARY_PERIODIC) || (by && by != DM_BOUNDARY_PERIODIC))) { 66347c6ae99SBarry Smith /* 66447c6ae99SBarry Smith Recompute the local to global mappings, this time keeping the 665ce00eea3SSatish Balay information about the cross corner processor numbers and any ghosted 666ce00eea3SSatish Balay but not periodic indices. 66747c6ae99SBarry Smith */ 66847c6ae99SBarry Smith nn = 0; 66947c6ae99SBarry Smith xbase = bases[rank]; 67047c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 67147c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 672ce00eea3SSatish Balay x_t = lx[n0 % m]; 67347c6ae99SBarry Smith y_t = ly[(n0/m)]; 67447c6ae99SBarry Smith s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x; 6758865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 676ce00eea3SSatish Balay } else if (xs-Xs > 0 && ys-Ys > 0) { 6778865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 67847c6ae99SBarry Smith } 67947c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 68047c6ae99SBarry Smith x_t = x; 68147c6ae99SBarry Smith y_t = ly[(n1/m)]; 68247c6ae99SBarry Smith s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t; 6838865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 684ce00eea3SSatish Balay } else if (ys-Ys > 0) { 685bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_MIRROR) { 6868865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(s_y - i + 1) + j; 687624904c4SBarry Smith } else { 6888865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 68947c6ae99SBarry Smith } 690624904c4SBarry Smith } 69147c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 692ce00eea3SSatish Balay x_t = lx[n2 % m]; 69347c6ae99SBarry Smith y_t = ly[(n2/m)]; 69447c6ae99SBarry Smith s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t; 6958865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 696ce00eea3SSatish Balay } else if (Xe-xe> 0 && ys-Ys > 0) { 6978865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 69847c6ae99SBarry Smith } 69947c6ae99SBarry Smith } 70047c6ae99SBarry Smith 70147c6ae99SBarry Smith for (i=0; i<y; i++) { 70247c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 703ce00eea3SSatish Balay x_t = lx[n3 % m]; 70447c6ae99SBarry Smith /* y_t = y; */ 70547c6ae99SBarry Smith s_t = bases[n3] + (i+1)*x_t - s_x; 7068865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 707ce00eea3SSatish Balay } else if (xs-Xs > 0) { 708bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_MIRROR) { 7098865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*i + s_x - j; 710624904c4SBarry Smith } else { 7118865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 71247c6ae99SBarry Smith } 713624904c4SBarry Smith } 71447c6ae99SBarry Smith 7158865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = xbase++; /* interior */ 71647c6ae99SBarry Smith 71747c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 718ce00eea3SSatish Balay x_t = lx[n5 % m]; 71947c6ae99SBarry Smith /* y_t = y; */ 72047c6ae99SBarry Smith s_t = bases[n5] + (i)*x_t; 7218865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 722ce00eea3SSatish Balay } else if (Xe-xe > 0) { 723bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_MIRROR) { 7248865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*(i + 1) - 2 - j; 725624904c4SBarry Smith } else { 7268865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 72747c6ae99SBarry Smith } 72847c6ae99SBarry Smith } 729624904c4SBarry Smith } 73047c6ae99SBarry Smith 73147c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 73247c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 733ce00eea3SSatish Balay x_t = lx[n6 % m]; 73447c6ae99SBarry Smith /* y_t = ly[(n6/m)]; */ 73547c6ae99SBarry Smith s_t = bases[n6] + (i)*x_t - s_x; 7368865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 737ce00eea3SSatish Balay } else if (xs-Xs > 0 && Ye-ye > 0) { 7388865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 73947c6ae99SBarry Smith } 74047c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 74147c6ae99SBarry Smith x_t = x; 74247c6ae99SBarry Smith /* y_t = ly[(n7/m)]; */ 74347c6ae99SBarry Smith s_t = bases[n7] + (i-1)*x_t; 7448865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 745ce00eea3SSatish Balay } else if (Ye-ye > 0) { 746bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_MIRROR) { 7478865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(y - i - 1) + j; 748624904c4SBarry Smith } else { 7498865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 75047c6ae99SBarry Smith } 751624904c4SBarry Smith } 75247c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 753ce00eea3SSatish Balay x_t = lx[n8 % m]; 75447c6ae99SBarry Smith /* y_t = ly[(n8/m)]; */ 75547c6ae99SBarry Smith s_t = bases[n8] + (i-1)*x_t; 7568865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 757ce00eea3SSatish Balay } else if (Xe-xe > 0 && Ye-ye > 0) { 7588865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 75947c6ae99SBarry Smith } 76047c6ae99SBarry Smith } 76147c6ae99SBarry Smith } 762ce00eea3SSatish Balay /* 763ce00eea3SSatish Balay Set the local to global ordering in the global vector, this allows use 764ce00eea3SSatish Balay of VecSetValuesLocal(). 765ce00eea3SSatish Balay */ 76645b6f7e9SBarry Smith ierr = ISLocalToGlobalMappingCreate(comm,dof,nn,idx,PETSC_OWN_POINTER,&da->ltogmap);CHKERRQ(ierr); 7673bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)da->ltogmap);CHKERRQ(ierr); 76847c6ae99SBarry Smith 769ce00eea3SSatish Balay ierr = PetscFree2(bases,ldims);CHKERRQ(ierr); 77047c6ae99SBarry Smith dd->m = m; dd->n = n; 771ce00eea3SSatish Balay /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */ 772ce00eea3SSatish Balay dd->xs = xs*dof; dd->xe = xe*dof; dd->ys = ys; dd->ye = ye; dd->zs = 0; dd->ze = 1; 773ce00eea3SSatish Balay dd->Xs = Xs*dof; dd->Xe = Xe*dof; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = 0; dd->Ze = 1; 77447c6ae99SBarry Smith 775fcfd50ebSBarry Smith ierr = VecDestroy(&local);CHKERRQ(ierr); 776fcfd50ebSBarry Smith ierr = VecDestroy(&global);CHKERRQ(ierr); 77747c6ae99SBarry Smith 77847c6ae99SBarry Smith dd->gtol = gtol; 77947c6ae99SBarry Smith dd->base = base; 7809a42bb27SBarry Smith da->ops->view = DMView_DA_2d; 7810298fd71SBarry Smith dd->ltol = NULL; 7820298fd71SBarry Smith dd->ao = NULL; 78347c6ae99SBarry Smith PetscFunctionReturn(0); 78447c6ae99SBarry Smith } 78547c6ae99SBarry Smith 78647c6ae99SBarry Smith /*@C 787aa219208SBarry Smith DMDACreate2d - Creates an object that will manage the communication of two-dimensional 78847c6ae99SBarry Smith regular array data that is distributed across some processors. 78947c6ae99SBarry Smith 79047c6ae99SBarry Smith Collective on MPI_Comm 79147c6ae99SBarry Smith 79247c6ae99SBarry Smith Input Parameters: 79347c6ae99SBarry Smith + comm - MPI communicator 7941321219cSEthan Coon . bx,by - type of ghost nodes the array have. 795bff4a2f0SMatthew G. Knepley Use one of DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_PERIODIC. 796aa219208SBarry Smith . stencil_type - stencil type. Use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR. 797897f7067SBarry Smith . M,N - global dimension in each direction of the array 79847c6ae99SBarry Smith . m,n - corresponding number of processors in each dimension 79947c6ae99SBarry Smith (or PETSC_DECIDE to have calculated) 80047c6ae99SBarry Smith . dof - number of degrees of freedom per node 80147c6ae99SBarry Smith . s - stencil width 80247c6ae99SBarry Smith - lx, ly - arrays containing the number of nodes in each cell along 8030298fd71SBarry Smith the x and y coordinates, or NULL. If non-null, these 80447c6ae99SBarry Smith must be of length as m and n, and the corresponding 80547c6ae99SBarry Smith m and n cannot be PETSC_DECIDE. The sum of the lx[] entries 80647c6ae99SBarry Smith must be M, and the sum of the ly[] entries must be N. 80747c6ae99SBarry Smith 80847c6ae99SBarry Smith Output Parameter: 80947c6ae99SBarry Smith . da - the resulting distributed array object 81047c6ae99SBarry Smith 81147c6ae99SBarry Smith Options Database Key: 812706a11cbSBarry Smith + -dm_view - Calls DMView() at the conclusion of DMDACreate2d() 81347c6ae99SBarry Smith . -da_grid_x <nx> - number of grid points in x direction, if M < 0 81447c6ae99SBarry Smith . -da_grid_y <ny> - number of grid points in y direction, if N < 0 81547c6ae99SBarry Smith . -da_processors_x <nx> - number of processors in x direction 81647c6ae99SBarry Smith . -da_processors_y <ny> - number of processors in y direction 817e0f5d30fSBarry Smith . -da_refine_x <rx> - refinement ratio in x direction 818e0f5d30fSBarry Smith . -da_refine_y <ry> - refinement ratio in y direction 819e0f5d30fSBarry Smith - -da_refine <n> - refine the DMDA n times before creating, if M or N < 0 820e0f5d30fSBarry Smith 82147c6ae99SBarry Smith 82247c6ae99SBarry Smith Level: beginner 82347c6ae99SBarry Smith 82447c6ae99SBarry Smith Notes: 825aa219208SBarry Smith The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the 826aa219208SBarry Smith standard 5-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes 82747c6ae99SBarry Smith the standard 9-pt stencil. 82847c6ae99SBarry Smith 829aa219208SBarry Smith The array data itself is NOT stored in the DMDA, it is stored in Vec objects; 830564755cdSBarry Smith The appropriate vector objects can be obtained with calls to DMCreateGlobalVector() 831564755cdSBarry Smith and DMCreateLocalVector() and calls to VecDuplicate() if more are needed. 83247c6ae99SBarry Smith 833897f7067SBarry Smith You must call DMSetUp() after this call before using this DM. 834897f7067SBarry Smith 835897f7067SBarry Smith If you wish to use the options database to change values in the DMDA call DMSetFromOptions() after this call 836897f7067SBarry Smith but before DMSetUp(). 837897f7067SBarry Smith 83847c6ae99SBarry Smith .keywords: distributed array, create, two-dimensional 83947c6ae99SBarry Smith 840aa219208SBarry Smith .seealso: DMDestroy(), DMView(), DMDACreate1d(), DMDACreate3d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(), 84199f0b487SRichard Tran Mills DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMLocalToLocalBegin(), DMLocalToLocalEnd(), DMDASetRefinementFactor(), 842d461ba97SBarry Smith DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges() 84347c6ae99SBarry Smith 84447c6ae99SBarry Smith @*/ 845fe16a2e9SBarry Smith 846bff4a2f0SMatthew G. Knepley PetscErrorCode DMDACreate2d(MPI_Comm comm,DMBoundaryType bx,DMBoundaryType by,DMDAStencilType stencil_type, 8479a42bb27SBarry Smith PetscInt M,PetscInt N,PetscInt m,PetscInt n,PetscInt dof,PetscInt s,const PetscInt lx[],const PetscInt ly[],DM *da) 84847c6ae99SBarry Smith { 84947c6ae99SBarry Smith PetscErrorCode ierr; 85047c6ae99SBarry Smith 85147c6ae99SBarry Smith PetscFunctionBegin; 852aa219208SBarry Smith ierr = DMDACreate(comm, da);CHKERRQ(ierr); 853c73cfb54SMatthew G. Knepley ierr = DMSetDimension(*da, 2);CHKERRQ(ierr); 854aa219208SBarry Smith ierr = DMDASetSizes(*da, M, N, 1);CHKERRQ(ierr); 855aa219208SBarry Smith ierr = DMDASetNumProcs(*da, m, n, PETSC_DECIDE);CHKERRQ(ierr); 856bff4a2f0SMatthew G. Knepley ierr = DMDASetBoundaryType(*da, bx, by, DM_BOUNDARY_NONE);CHKERRQ(ierr); 857aa219208SBarry Smith ierr = DMDASetDof(*da, dof);CHKERRQ(ierr); 858aa219208SBarry Smith ierr = DMDASetStencilType(*da, stencil_type);CHKERRQ(ierr); 859aa219208SBarry Smith ierr = DMDASetStencilWidth(*da, s);CHKERRQ(ierr); 8600298fd71SBarry Smith ierr = DMDASetOwnershipRanges(*da, lx, ly, NULL);CHKERRQ(ierr); 86147c6ae99SBarry Smith PetscFunctionReturn(0); 86247c6ae99SBarry Smith } 863