19a42bb27SBarry Smith 2b45d2f2cSJed Brown #include <petsc-private/daimpl.h> /*I "petscdmda.h" I*/ 347c6ae99SBarry Smith 447c6ae99SBarry Smith #undef __FUNCT__ 59a42bb27SBarry Smith #define __FUNCT__ "DMView_DA_2d" 69a42bb27SBarry Smith PetscErrorCode DMView_DA_2d(DM da,PetscViewer viewer) 747c6ae99SBarry Smith { 847c6ae99SBarry Smith PetscErrorCode ierr; 947c6ae99SBarry Smith PetscMPIInt rank; 109a42bb27SBarry Smith PetscBool iascii,isdraw,isbinary; 1147c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 129a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 139a42bb27SBarry Smith PetscBool ismatlab; 149a42bb27SBarry Smith #endif 1547c6ae99SBarry Smith 1647c6ae99SBarry Smith PetscFunctionBegin; 1747c6ae99SBarry Smith ierr = MPI_Comm_rank(((PetscObject)da)->comm,&rank);CHKERRQ(ierr); 1847c6ae99SBarry Smith 19251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 20251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);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); 2947c6ae99SBarry Smith if (format != PETSC_VIEWER_ASCII_VTK && format != PETSC_VIEWER_ASCII_VTK_CELL) { 30aa219208SBarry Smith DMDALocalInfo info; 31aa219208SBarry Smith ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr); 327b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);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); 367b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr); 373da9ae13SJed Brown } else { 383da9ae13SJed Brown ierr = DMView_DA_VTK(da,viewer);CHKERRQ(ierr); 3947c6ae99SBarry Smith } 4047c6ae99SBarry Smith } else if (isdraw) { 4147c6ae99SBarry Smith PetscDraw draw; 4247c6ae99SBarry Smith double ymin = -1*dd->s-1,ymax = dd->N+dd->s; 4347c6ae99SBarry Smith double xmin = -1*dd->s-1,xmax = dd->M+dd->s; 4447c6ae99SBarry Smith double x,y; 4547c6ae99SBarry Smith PetscInt base,*idx; 4647c6ae99SBarry Smith char node[10]; 4747c6ae99SBarry Smith PetscBool isnull; 4847c6ae99SBarry Smith 4947c6ae99SBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 5047c6ae99SBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 516636e97aSMatthew G Knepley if (!da->coordinates) { 5247c6ae99SBarry Smith ierr = PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);CHKERRQ(ierr); 5347c6ae99SBarry Smith } 5447c6ae99SBarry Smith ierr = PetscDrawSynchronizedClear(draw);CHKERRQ(ierr); 5547c6ae99SBarry Smith 5647c6ae99SBarry Smith /* first processor draw all node lines */ 5747c6ae99SBarry Smith if (!rank) { 5847c6ae99SBarry Smith ymin = 0.0; ymax = dd->N - 1; 5947c6ae99SBarry Smith for (xmin=0; xmin<dd->M; xmin++) { 6047c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr); 6147c6ae99SBarry Smith } 6247c6ae99SBarry Smith xmin = 0.0; xmax = dd->M - 1; 6347c6ae99SBarry Smith for (ymin=0; ymin<dd->N; ymin++) { 6447c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);CHKERRQ(ierr); 6547c6ae99SBarry Smith } 6647c6ae99SBarry Smith } 6747c6ae99SBarry Smith ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr); 6847c6ae99SBarry Smith ierr = PetscDrawPause(draw);CHKERRQ(ierr); 6947c6ae99SBarry Smith 7047c6ae99SBarry Smith /* draw my box */ 7147c6ae99SBarry Smith ymin = dd->ys; ymax = dd->ye - 1; xmin = dd->xs/dd->w; 7247c6ae99SBarry Smith xmax =(dd->xe-1)/dd->w; 7347c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);CHKERRQ(ierr); 7447c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 7547c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 7647c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 7747c6ae99SBarry Smith 7847c6ae99SBarry Smith /* put in numbers */ 7947c6ae99SBarry Smith base = (dd->base)/dd->w; 8047c6ae99SBarry Smith for (y=ymin; y<=ymax; y++) { 8147c6ae99SBarry Smith for (x=xmin; x<=xmax; x++) { 8247c6ae99SBarry Smith sprintf(node,"%d",(int)base++); 8347c6ae99SBarry Smith ierr = PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node);CHKERRQ(ierr); 8447c6ae99SBarry Smith } 8547c6ae99SBarry Smith } 8647c6ae99SBarry Smith 8747c6ae99SBarry Smith ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr); 8847c6ae99SBarry Smith ierr = PetscDrawPause(draw);CHKERRQ(ierr); 8947c6ae99SBarry Smith /* overlay ghost numbers, useful for error checking */ 9047c6ae99SBarry Smith /* put in numbers */ 9147c6ae99SBarry Smith 9247c6ae99SBarry Smith base = 0; idx = dd->idx; 9347c6ae99SBarry Smith ymin = dd->Ys; ymax = dd->Ye; xmin = dd->Xs; xmax = dd->Xe; 9447c6ae99SBarry Smith for (y=ymin; y<ymax; y++) { 9547c6ae99SBarry Smith for (x=xmin; x<xmax; x++) { 9647c6ae99SBarry Smith if ((base % dd->w) == 0) { 9747c6ae99SBarry Smith sprintf(node,"%d",(int)(idx[base]/dd->w)); 9847c6ae99SBarry Smith ierr = PetscDrawString(draw,x/dd->w,y,PETSC_DRAW_BLUE,node);CHKERRQ(ierr); 9947c6ae99SBarry Smith } 10047c6ae99SBarry Smith base++; 10147c6ae99SBarry Smith } 10247c6ae99SBarry Smith } 10347c6ae99SBarry Smith ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr); 10447c6ae99SBarry Smith ierr = PetscDrawPause(draw);CHKERRQ(ierr); 1059a42bb27SBarry Smith } else if (isbinary){ 1069a42bb27SBarry Smith ierr = DMView_DA_Binary(da,viewer);CHKERRQ(ierr); 1079a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 1089a42bb27SBarry Smith } else if (ismatlab) { 1099a42bb27SBarry Smith ierr = DMView_DA_Matlab(da,viewer);CHKERRQ(ierr); 1109a42bb27SBarry Smith #endif 111aa219208SBarry Smith } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for DMDA 1d",((PetscObject)viewer)->type_name); 11247c6ae99SBarry Smith PetscFunctionReturn(0); 11347c6ae99SBarry Smith } 11447c6ae99SBarry Smith 11547c6ae99SBarry Smith /* 11647c6ae99SBarry Smith M is number of grid points 11747c6ae99SBarry Smith m is number of processors 11847c6ae99SBarry Smith 11947c6ae99SBarry Smith */ 12047c6ae99SBarry Smith #undef __FUNCT__ 121aa219208SBarry Smith #define __FUNCT__ "DMDASplitComm2d" 1227087cfbeSBarry Smith PetscErrorCode DMDASplitComm2d(MPI_Comm comm,PetscInt M,PetscInt N,PetscInt sw,MPI_Comm *outcomm) 12347c6ae99SBarry Smith { 12447c6ae99SBarry Smith PetscErrorCode ierr; 12547c6ae99SBarry Smith PetscInt m,n = 0,x = 0,y = 0; 12647c6ae99SBarry Smith PetscMPIInt size,csize,rank; 12747c6ae99SBarry Smith 12847c6ae99SBarry Smith PetscFunctionBegin; 12947c6ae99SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 13047c6ae99SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 13147c6ae99SBarry Smith 13247c6ae99SBarry Smith csize = 4*size; 13347c6ae99SBarry Smith do { 13447c6ae99SBarry 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); 13547c6ae99SBarry Smith csize = csize/4; 13647c6ae99SBarry Smith 13747c6ae99SBarry Smith m = (PetscInt)(0.5 + sqrt(((double)M)*((double)csize)/((double)N))); 13847c6ae99SBarry Smith if (!m) m = 1; 13947c6ae99SBarry Smith while (m > 0) { 14047c6ae99SBarry Smith n = csize/m; 14147c6ae99SBarry Smith if (m*n == csize) break; 14247c6ae99SBarry Smith m--; 14347c6ae99SBarry Smith } 14447c6ae99SBarry Smith if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;} 14547c6ae99SBarry Smith 14647c6ae99SBarry Smith x = M/m + ((M % m) > ((csize-1) % m)); 14747c6ae99SBarry Smith y = (N + (csize-1)/m)/n; 14847c6ae99SBarry Smith } while ((x < 4 || y < 4) && csize > 1); 14947c6ae99SBarry Smith if (size != csize) { 15047c6ae99SBarry Smith MPI_Group entire_group,sub_group; 15147c6ae99SBarry Smith PetscMPIInt i,*groupies; 15247c6ae99SBarry Smith 15347c6ae99SBarry Smith ierr = MPI_Comm_group(comm,&entire_group);CHKERRQ(ierr); 15447c6ae99SBarry Smith ierr = PetscMalloc(csize*sizeof(PetscInt),&groupies);CHKERRQ(ierr); 15547c6ae99SBarry Smith for (i=0; i<csize; i++) { 15647c6ae99SBarry Smith groupies[i] = (rank/csize)*csize + i; 15747c6ae99SBarry Smith } 15847c6ae99SBarry Smith ierr = MPI_Group_incl(entire_group,csize,groupies,&sub_group);CHKERRQ(ierr); 15947c6ae99SBarry Smith ierr = PetscFree(groupies);CHKERRQ(ierr); 16047c6ae99SBarry Smith ierr = MPI_Comm_create(comm,sub_group,outcomm);CHKERRQ(ierr); 16147c6ae99SBarry Smith ierr = MPI_Group_free(&entire_group);CHKERRQ(ierr); 16247c6ae99SBarry Smith ierr = MPI_Group_free(&sub_group);CHKERRQ(ierr); 163aa219208SBarry Smith ierr = PetscInfo1(0,"DMDASplitComm2d:Creating redundant coarse problems of size %d\n",csize);CHKERRQ(ierr); 16447c6ae99SBarry Smith } else { 16547c6ae99SBarry Smith *outcomm = comm; 16647c6ae99SBarry Smith } 16747c6ae99SBarry Smith PetscFunctionReturn(0); 16847c6ae99SBarry Smith } 16947c6ae99SBarry Smith 17047c6ae99SBarry Smith #if defined(new) 17147c6ae99SBarry Smith #undef __FUNCT__ 172aa219208SBarry Smith #define __FUNCT__ "DMDAGetDiagonal_MFFD" 17347c6ae99SBarry Smith /* 174aa219208SBarry Smith DMDAGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix where local 175aa219208SBarry Smith function lives on a DMDA 17647c6ae99SBarry Smith 17747c6ae99SBarry Smith y ~= (F(u + ha) - F(u))/h, 17847c6ae99SBarry Smith where F = nonlinear function, as set by SNESSetFunction() 17947c6ae99SBarry Smith u = current iterate 18047c6ae99SBarry Smith h = difference interval 18147c6ae99SBarry Smith */ 182aa219208SBarry Smith PetscErrorCode DMDAGetDiagonal_MFFD(DM da,Vec U,Vec a) 18347c6ae99SBarry Smith { 18447c6ae99SBarry Smith PetscScalar h,*aa,*ww,v; 18547c6ae99SBarry Smith PetscReal epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON; 18647c6ae99SBarry Smith PetscErrorCode ierr; 18747c6ae99SBarry Smith PetscInt gI,nI; 18847c6ae99SBarry Smith MatStencil stencil; 189aa219208SBarry Smith DMDALocalInfo info; 19047c6ae99SBarry Smith 19147c6ae99SBarry Smith PetscFunctionBegin; 19247c6ae99SBarry Smith ierr = (*ctx->func)(0,U,a,ctx->funcctx);CHKERRQ(ierr); 19347c6ae99SBarry Smith ierr = (*ctx->funcisetbase)(U,ctx->funcctx);CHKERRQ(ierr); 19447c6ae99SBarry Smith 19547c6ae99SBarry Smith ierr = VecGetArray(U,&ww);CHKERRQ(ierr); 19647c6ae99SBarry Smith ierr = VecGetArray(a,&aa);CHKERRQ(ierr); 19747c6ae99SBarry Smith 19847c6ae99SBarry Smith nI = 0; 19947c6ae99SBarry Smith h = ww[gI]; 20047c6ae99SBarry Smith if (h == 0.0) h = 1.0; 20147c6ae99SBarry Smith if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin; 20247c6ae99SBarry Smith else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin; 20347c6ae99SBarry Smith h *= epsilon; 20447c6ae99SBarry Smith 20547c6ae99SBarry Smith ww[gI] += h; 20647c6ae99SBarry Smith ierr = (*ctx->funci)(i,w,&v,ctx->funcctx);CHKERRQ(ierr); 20747c6ae99SBarry Smith aa[nI] = (v - aa[nI])/h; 20847c6ae99SBarry Smith ww[gI] -= h; 20947c6ae99SBarry Smith nI++; 21047c6ae99SBarry Smith } 21147c6ae99SBarry Smith ierr = VecRestoreArray(U,&ww);CHKERRQ(ierr); 21247c6ae99SBarry Smith ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr); 21347c6ae99SBarry Smith PetscFunctionReturn(0); 21447c6ae99SBarry Smith } 21547c6ae99SBarry Smith #endif 21647c6ae99SBarry Smith 21747c6ae99SBarry Smith #undef __FUNCT__ 2189a42bb27SBarry Smith #define __FUNCT__ "DMSetUp_DA_2D" 2197087cfbeSBarry Smith PetscErrorCode DMSetUp_DA_2D(DM da) 22047c6ae99SBarry Smith { 22147c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 22247c6ae99SBarry Smith const PetscInt M = dd->M; 22347c6ae99SBarry Smith const PetscInt N = dd->N; 22447c6ae99SBarry Smith PetscInt m = dd->m; 22547c6ae99SBarry Smith PetscInt n = dd->n; 22688661749SPeter Brune PetscInt o = dd->overlap; 22747c6ae99SBarry Smith const PetscInt dof = dd->w; 22847c6ae99SBarry Smith const PetscInt s = dd->s; 22919fd82e9SBarry Smith DMDABoundaryType bx = dd->bx; 23019fd82e9SBarry Smith DMDABoundaryType by = dd->by; 23119fd82e9SBarry Smith DMDAStencilType stencil_type = dd->stencil_type; 23247c6ae99SBarry Smith PetscInt *lx = dd->lx; 23347c6ae99SBarry Smith PetscInt *ly = dd->ly; 23447c6ae99SBarry Smith MPI_Comm comm; 23547c6ae99SBarry Smith PetscMPIInt rank,size; 236ce00eea3SSatish Balay PetscInt xs,xe,ys,ye,x,y,Xs,Xe,Ys,Ye,start,end,IXs,IXe,IYs,IYe; 237ce00eea3SSatish Balay PetscInt up,down,left,right,i,n0,n1,n2,n3,n5,n6,n7,n8,*idx,nn,*idx_cpy; 238ce00eea3SSatish Balay const PetscInt *idx_full; 239db87c5ecSEthan Coon PetscInt xbase,*bases,*ldims,j,x_t,y_t,s_t,base,count; 24047c6ae99SBarry Smith PetscInt s_x,s_y; /* s proportionalized to w */ 24147c6ae99SBarry Smith PetscInt sn0 = 0,sn2 = 0,sn6 = 0,sn8 = 0; 24247c6ae99SBarry Smith Vec local,global; 24347c6ae99SBarry Smith VecScatter ltog,gtol; 244ce00eea3SSatish Balay IS to,from,ltogis; 24547c6ae99SBarry Smith PetscErrorCode ierr; 24647c6ae99SBarry Smith 24747c6ae99SBarry Smith PetscFunctionBegin; 248cc76e2cfSBarry Smith if (stencil_type == DMDA_STENCIL_BOX && (bx == DMDA_BOUNDARY_MIRROR || by == DMDA_BOUNDARY_MIRROR)) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Mirror boundary and box stencil"); 24947c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 2503855c12bSBarry Smith #if !defined(PETSC_USE_64BIT_INDICES) 2513855c12bSBarry 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); 2523855c12bSBarry Smith #endif 25347c6ae99SBarry Smith 25447c6ae99SBarry Smith if (dof < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have 1 or more degrees of freedom per node: %D",dof); 25547c6ae99SBarry Smith if (s < 0) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Stencil width cannot be negative: %D",s); 25647c6ae99SBarry Smith 25747c6ae99SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 25847c6ae99SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 25947c6ae99SBarry Smith 26047c6ae99SBarry Smith if (m != PETSC_DECIDE) { 26147c6ae99SBarry Smith if (m < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m); 26247c6ae99SBarry Smith else if (m > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size); 26347c6ae99SBarry Smith } 26447c6ae99SBarry Smith if (n != PETSC_DECIDE) { 26547c6ae99SBarry Smith if (n < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n); 26647c6ae99SBarry Smith else if (n > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size); 26747c6ae99SBarry Smith } 26847c6ae99SBarry Smith 26947c6ae99SBarry Smith if (m == PETSC_DECIDE || n == PETSC_DECIDE) { 27047c6ae99SBarry Smith if (n != PETSC_DECIDE) { 27147c6ae99SBarry Smith m = size/n; 27247c6ae99SBarry Smith } else if (m != PETSC_DECIDE) { 27347c6ae99SBarry Smith n = size/m; 27447c6ae99SBarry Smith } else { 27547c6ae99SBarry Smith /* try for squarish distribution */ 27647c6ae99SBarry Smith m = (PetscInt)(0.5 + sqrt(((double)M)*((double)size)/((double)N))); 27747c6ae99SBarry Smith if (!m) m = 1; 27847c6ae99SBarry Smith while (m > 0) { 27947c6ae99SBarry Smith n = size/m; 28047c6ae99SBarry Smith if (m*n == size) break; 28147c6ae99SBarry Smith m--; 28247c6ae99SBarry Smith } 28347c6ae99SBarry Smith if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;} 28447c6ae99SBarry Smith } 28547c6ae99SBarry 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 "); 28647c6ae99SBarry Smith } else if (m*n != size) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition"); 28747c6ae99SBarry Smith 28847c6ae99SBarry Smith if (M < m) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m); 28947c6ae99SBarry Smith if (N < n) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n); 29047c6ae99SBarry Smith 29147c6ae99SBarry Smith /* 29247c6ae99SBarry Smith Determine locally owned region 29347c6ae99SBarry Smith xs is the first local node number, x is the number of local nodes 29447c6ae99SBarry Smith */ 29547c6ae99SBarry Smith if (!lx) { 29647c6ae99SBarry Smith ierr = PetscMalloc(m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr); 29747c6ae99SBarry Smith lx = dd->lx; 29847c6ae99SBarry Smith for (i=0; i<m; i++) { 29947c6ae99SBarry Smith lx[i] = M/m + ((M % m) > i); 30047c6ae99SBarry Smith } 30147c6ae99SBarry Smith } 30247c6ae99SBarry Smith x = lx[rank % m]; 30347c6ae99SBarry Smith xs = 0; 30447c6ae99SBarry Smith for (i=0; i<(rank % m); i++) { 30547c6ae99SBarry Smith xs += lx[i]; 30647c6ae99SBarry Smith } 30747c6ae99SBarry Smith #if defined(PETSC_USE_DEBUG) 30847c6ae99SBarry Smith left = xs; 30947c6ae99SBarry Smith for (i=(rank % m); i<m; i++) { 31047c6ae99SBarry Smith left += lx[i]; 31147c6ae99SBarry Smith } 31247c6ae99SBarry 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); 31347c6ae99SBarry Smith #endif 31447c6ae99SBarry Smith 31547c6ae99SBarry Smith /* 31647c6ae99SBarry Smith Determine locally owned region 31747c6ae99SBarry Smith ys is the first local node number, y is the number of local nodes 31847c6ae99SBarry Smith */ 31947c6ae99SBarry Smith if (!ly) { 32047c6ae99SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt), &dd->ly);CHKERRQ(ierr); 32147c6ae99SBarry Smith ly = dd->ly; 32247c6ae99SBarry Smith for (i=0; i<n; i++) { 32347c6ae99SBarry Smith ly[i] = N/n + ((N % n) > i); 32447c6ae99SBarry Smith } 32547c6ae99SBarry Smith } 32647c6ae99SBarry Smith y = ly[rank/m]; 32747c6ae99SBarry Smith ys = 0; 32847c6ae99SBarry Smith for (i=0; i<(rank/m); i++) { 32947c6ae99SBarry Smith ys += ly[i]; 33047c6ae99SBarry Smith } 33147c6ae99SBarry Smith #if defined(PETSC_USE_DEBUG) 33247c6ae99SBarry Smith left = ys; 33347c6ae99SBarry Smith for (i=(rank/m); i<n; i++) { 33447c6ae99SBarry Smith left += ly[i]; 33547c6ae99SBarry Smith } 33647c6ae99SBarry 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); 33747c6ae99SBarry Smith #endif 33847c6ae99SBarry Smith 339bcea557cSEthan Coon /* 340bcea557cSEthan Coon check if the scatter requires more than one process neighbor or wraps around 341bcea557cSEthan Coon the domain more than once 342bcea557cSEthan Coon */ 34388661749SPeter Brune if ((x < s+o) && ((m > 1) || (bx == DMDA_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+o); 34488661749SPeter Brune if ((y < s+o) && ((n > 1) || (by == DMDA_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+o); 34547c6ae99SBarry Smith xe = xs + x; 34647c6ae99SBarry Smith ye = ys + y; 34747c6ae99SBarry Smith 348ce00eea3SSatish Balay /* determine ghost region (Xs) and region scattered into (IXs) */ 34988661749SPeter Brune if (xs-s-o > 0) { 35088661749SPeter Brune Xs = xs - s - o; IXs = xs - s - o; 35188661749SPeter Brune } else { 35288661749SPeter Brune if (bx) { 35388661749SPeter Brune Xs = xs - s; 35488661749SPeter Brune } else { 35588661749SPeter Brune Xs = 0; 35688661749SPeter Brune } 35788661749SPeter Brune IXs = 0; 35888661749SPeter Brune } 35988661749SPeter Brune if (xe+s+o <= M) { 36088661749SPeter Brune Xe = xe + s + o; IXe = xe + s + o; 36188661749SPeter Brune } else { 36288661749SPeter Brune if (bx) { 36388661749SPeter Brune Xs = xs - s - o; Xe = xe + s; 36488661749SPeter Brune } else { 36588661749SPeter Brune Xe = M; 36688661749SPeter Brune } 36788661749SPeter Brune IXe = M; 36888661749SPeter Brune } 36947c6ae99SBarry Smith 370c2859e5eSBarry Smith if (bx == DMDA_BOUNDARY_PERIODIC || bx == DMDA_BOUNDARY_MIRROR) { 37188661749SPeter Brune IXs = xs - s - o; 37288661749SPeter Brune IXe = xe + s + o; 37388661749SPeter Brune Xs = xs - s - o; 37488661749SPeter Brune Xe = xe + s + o; 37588661749SPeter Brune } 37647c6ae99SBarry Smith 37788661749SPeter Brune if (ys-s-o > 0) { 37888661749SPeter Brune Ys = ys - s - o; IYs = ys - s - o; 37988661749SPeter Brune } else { 38088661749SPeter Brune if (by) { 38188661749SPeter Brune Ys = ys - s; 38288661749SPeter Brune } else { 38388661749SPeter Brune Ys = 0; 38488661749SPeter Brune } 38588661749SPeter Brune IYs = 0; 38688661749SPeter Brune } 38788661749SPeter Brune if (ye+s+o <= N) { 38888661749SPeter Brune Ye = ye + s + o; IYe = ye + s + o; 38988661749SPeter Brune } else { 39088661749SPeter Brune if (by) { 39188661749SPeter Brune Ye = ye + s; 39288661749SPeter Brune } else { 39388661749SPeter Brune Ye = N; 39488661749SPeter Brune } 39588661749SPeter Brune IYe = N; 39688661749SPeter Brune } 39788661749SPeter Brune 398ac119b13SBarry Smith if (by == DMDA_BOUNDARY_PERIODIC || by == DMDA_BOUNDARY_MIRROR) { 39988661749SPeter Brune IYs = ys - s - o; 40088661749SPeter Brune IYe = ye + s + o; 40188661749SPeter Brune Ys = ys - s - o; 40288661749SPeter Brune Ye = ye + s + o; 40388661749SPeter Brune } 40488661749SPeter Brune 40588661749SPeter Brune /* stencil length in each direction */ 40688661749SPeter Brune s_x = s + o; 40788661749SPeter Brune s_y = s + o; 40847c6ae99SBarry Smith 40947c6ae99SBarry Smith /* determine starting point of each processor */ 41047c6ae99SBarry Smith nn = x*y; 41147c6ae99SBarry Smith ierr = PetscMalloc2(size+1,PetscInt,&bases,size,PetscInt,&ldims);CHKERRQ(ierr); 41247c6ae99SBarry Smith ierr = MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);CHKERRQ(ierr); 41347c6ae99SBarry Smith bases[0] = 0; 41447c6ae99SBarry Smith for (i=1; i<=size; i++) { 41547c6ae99SBarry Smith bases[i] = ldims[i-1]; 41647c6ae99SBarry Smith } 41747c6ae99SBarry Smith for (i=1; i<=size; i++) { 41847c6ae99SBarry Smith bases[i] += bases[i-1]; 41947c6ae99SBarry Smith } 420ce00eea3SSatish Balay base = bases[rank]*dof; 42147c6ae99SBarry Smith 42247c6ae99SBarry Smith /* allocate the base parallel and sequential vectors */ 423ce00eea3SSatish Balay dd->Nlocal = x*y*dof; 424778a2246SBarry Smith ierr = VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,0,&global);CHKERRQ(ierr); 425ce00eea3SSatish Balay dd->nlocal = (Xe-Xs)*(Ye-Ys)*dof; 426778a2246SBarry Smith ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,0,&local);CHKERRQ(ierr); 42747c6ae99SBarry Smith 42847c6ae99SBarry Smith /* generate appropriate vector scatters */ 42947c6ae99SBarry Smith /* local to global inserts non-ghost point region into global */ 43047c6ae99SBarry Smith ierr = VecGetOwnershipRange(global,&start,&end);CHKERRQ(ierr); 431ce00eea3SSatish Balay ierr = ISCreateStride(comm,x*y*dof,start,1,&to);CHKERRQ(ierr); 43247c6ae99SBarry Smith 433ce00eea3SSatish Balay ierr = PetscMalloc(x*y*sizeof(PetscInt),&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 ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&from);CHKERRQ(ierr); 44447c6ae99SBarry Smith ierr = VecScatterCreate(local,from,global,to,<og);CHKERRQ(ierr); 44547c6ae99SBarry Smith ierr = PetscLogObjectParent(dd,ltog);CHKERRQ(ierr); 446fcfd50ebSBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 447fcfd50ebSBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 44847c6ae99SBarry Smith 449ce00eea3SSatish Balay /* global to local must include ghost points within the domain, 450ce00eea3SSatish Balay but not ghost points outside the domain that aren't periodic */ 45188661749SPeter Brune if (stencil_type == DMDA_STENCIL_BOX || o > 0) { 452db87c5ecSEthan Coon count = (IXe-IXs)*(IYe-IYs); 453db87c5ecSEthan Coon ierr = PetscMalloc(count*sizeof(PetscInt),&idx);CHKERRQ(ierr); 454ce00eea3SSatish Balay 455ce00eea3SSatish Balay left = IXs - Xs; right = left + (IXe-IXs); 456ce00eea3SSatish Balay down = IYs - Ys; up = down + (IYe-IYs); 457ce00eea3SSatish Balay count = 0; 458ce00eea3SSatish Balay for (i=down; i<up; i++) { 459ce00eea3SSatish Balay for (j=left; j<right; j++) { 460ce00eea3SSatish Balay idx[count++] = j + i*(Xe-Xs); 461ce00eea3SSatish Balay } 462ce00eea3SSatish Balay } 463ce00eea3SSatish Balay ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr); 464ce00eea3SSatish Balay 46547c6ae99SBarry Smith } else { 46647c6ae99SBarry Smith /* must drop into cross shape region */ 46747c6ae99SBarry Smith /* ---------| 46847c6ae99SBarry Smith | top | 469ce00eea3SSatish Balay |--- ---| up 47047c6ae99SBarry Smith | middle | 47147c6ae99SBarry Smith | | 472ce00eea3SSatish Balay ---- ---- down 47347c6ae99SBarry Smith | bottom | 47447c6ae99SBarry Smith ----------- 47547c6ae99SBarry Smith Xs xs xe Xe */ 476db87c5ecSEthan Coon count = (ys-IYs)*x + y*(IXe-IXs) + (IYe-ye)*x; 477db87c5ecSEthan Coon ierr = PetscMalloc(count*sizeof(PetscInt),&idx);CHKERRQ(ierr); 478ce00eea3SSatish Balay 479ce00eea3SSatish Balay left = xs - Xs; right = left + x; 480ce00eea3SSatish Balay down = ys - Ys; up = down + y; 48147c6ae99SBarry Smith count = 0; 482ce00eea3SSatish Balay /* bottom */ 483ce00eea3SSatish Balay for (i=(IYs-Ys); i<down; i++) { 484ce00eea3SSatish Balay for (j=left; j<right; j++) { 485ce00eea3SSatish Balay idx[count++] = j + i*(Xe-Xs); 48647c6ae99SBarry Smith } 48747c6ae99SBarry Smith } 48847c6ae99SBarry Smith /* middle */ 48947c6ae99SBarry Smith for (i=down; i<up; i++) { 490ce00eea3SSatish Balay for (j=(IXs-Xs); j<(IXe-Xs); j++) { 491ce00eea3SSatish Balay idx[count++] = j + i*(Xe-Xs); 49247c6ae99SBarry Smith } 49347c6ae99SBarry Smith } 49447c6ae99SBarry Smith /* top */ 495ce00eea3SSatish Balay for (i=up; i<up+IYe-ye; i++) { 496ce00eea3SSatish Balay for (j=left; j<right; j++) { 497ce00eea3SSatish Balay idx[count++] = j + i*(Xe-Xs); 49847c6ae99SBarry Smith } 49947c6ae99SBarry Smith } 50047c6ae99SBarry Smith ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr); 50147c6ae99SBarry Smith } 50247c6ae99SBarry Smith 50347c6ae99SBarry Smith 50447c6ae99SBarry Smith /* determine who lies on each side of us stored in n6 n7 n8 50547c6ae99SBarry Smith n3 n5 50647c6ae99SBarry Smith n0 n1 n2 50747c6ae99SBarry Smith */ 50847c6ae99SBarry Smith 50947c6ae99SBarry Smith /* Assume the Non-Periodic Case */ 51047c6ae99SBarry Smith n1 = rank - m; 51147c6ae99SBarry Smith if (rank % m) { 51247c6ae99SBarry Smith n0 = n1 - 1; 51347c6ae99SBarry Smith } else { 51447c6ae99SBarry Smith n0 = -1; 51547c6ae99SBarry Smith } 51647c6ae99SBarry Smith if ((rank+1) % m) { 51747c6ae99SBarry Smith n2 = n1 + 1; 51847c6ae99SBarry Smith n5 = rank + 1; 51947c6ae99SBarry Smith n8 = rank + m + 1; if (n8 >= m*n) n8 = -1; 52047c6ae99SBarry Smith } else { 52147c6ae99SBarry Smith n2 = -1; n5 = -1; n8 = -1; 52247c6ae99SBarry Smith } 52347c6ae99SBarry Smith if (rank % m) { 52447c6ae99SBarry Smith n3 = rank - 1; 52547c6ae99SBarry Smith n6 = n3 + m; if (n6 >= m*n) n6 = -1; 52647c6ae99SBarry Smith } else { 52747c6ae99SBarry Smith n3 = -1; n6 = -1; 52847c6ae99SBarry Smith } 52947c6ae99SBarry Smith n7 = rank + m; if (n7 >= m*n) n7 = -1; 53047c6ae99SBarry Smith 5311321219cSEthan Coon if (bx == DMDA_BOUNDARY_PERIODIC && by == DMDA_BOUNDARY_PERIODIC) { 53247c6ae99SBarry Smith /* Modify for Periodic Cases */ 53347c6ae99SBarry Smith /* Handle all four corners */ 53447c6ae99SBarry Smith if ((n6 < 0) && (n7 < 0) && (n3 < 0)) n6 = m-1; 53547c6ae99SBarry Smith if ((n8 < 0) && (n7 < 0) && (n5 < 0)) n8 = 0; 53647c6ae99SBarry Smith if ((n2 < 0) && (n5 < 0) && (n1 < 0)) n2 = size-m; 53747c6ae99SBarry Smith if ((n0 < 0) && (n3 < 0) && (n1 < 0)) n0 = size-1; 53847c6ae99SBarry Smith 53947c6ae99SBarry Smith /* Handle Top and Bottom Sides */ 54047c6ae99SBarry Smith if (n1 < 0) n1 = rank + m * (n-1); 54147c6ae99SBarry Smith if (n7 < 0) n7 = rank - m * (n-1); 54247c6ae99SBarry Smith if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1; 54347c6ae99SBarry Smith if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1; 54447c6ae99SBarry Smith if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1; 54547c6ae99SBarry Smith if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1; 54647c6ae99SBarry Smith 54747c6ae99SBarry Smith /* Handle Left and Right Sides */ 54847c6ae99SBarry Smith if (n3 < 0) n3 = rank + (m-1); 54947c6ae99SBarry Smith if (n5 < 0) n5 = rank - (m-1); 55047c6ae99SBarry Smith if ((n1 >= 0) && (n0 < 0)) n0 = rank-1; 55147c6ae99SBarry Smith if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1; 55247c6ae99SBarry Smith if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1; 55347c6ae99SBarry Smith if ((n7 >= 0) && (n8 < 0)) n8 = rank+1; 5541321219cSEthan Coon } else if (by == DMDA_BOUNDARY_PERIODIC) { /* Handle Top and Bottom Sides */ 555ce00eea3SSatish Balay if (n1 < 0) n1 = rank + m * (n-1); 556ce00eea3SSatish Balay if (n7 < 0) n7 = rank - m * (n-1); 557ce00eea3SSatish Balay if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1; 558ce00eea3SSatish Balay if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1; 559ce00eea3SSatish Balay if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1; 560ce00eea3SSatish Balay if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1; 5611321219cSEthan Coon } else if (bx == DMDA_BOUNDARY_PERIODIC) { /* Handle Left and Right Sides */ 562ce00eea3SSatish Balay if (n3 < 0) n3 = rank + (m-1); 563ce00eea3SSatish Balay if (n5 < 0) n5 = rank - (m-1); 564ce00eea3SSatish Balay if ((n1 >= 0) && (n0 < 0)) n0 = rank-1; 565ce00eea3SSatish Balay if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1; 566ce00eea3SSatish Balay if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1; 567ce00eea3SSatish Balay if ((n7 >= 0) && (n8 < 0)) n8 = rank+1; 56847c6ae99SBarry Smith } 569ce00eea3SSatish Balay 57047c6ae99SBarry Smith ierr = PetscMalloc(9*sizeof(PetscInt),&dd->neighbors);CHKERRQ(ierr); 57147c6ae99SBarry Smith dd->neighbors[0] = n0; 57247c6ae99SBarry Smith dd->neighbors[1] = n1; 57347c6ae99SBarry Smith dd->neighbors[2] = n2; 57447c6ae99SBarry Smith dd->neighbors[3] = n3; 57547c6ae99SBarry Smith dd->neighbors[4] = rank; 57647c6ae99SBarry Smith dd->neighbors[5] = n5; 57747c6ae99SBarry Smith dd->neighbors[6] = n6; 57847c6ae99SBarry Smith dd->neighbors[7] = n7; 57947c6ae99SBarry Smith dd->neighbors[8] = n8; 58047c6ae99SBarry Smith 58188661749SPeter Brune if (stencil_type == DMDA_STENCIL_STAR && o == 0) { 58247c6ae99SBarry Smith /* save corner processor numbers */ 58347c6ae99SBarry Smith sn0 = n0; sn2 = n2; sn6 = n6; sn8 = n8; 58447c6ae99SBarry Smith n0 = n2 = n6 = n8 = -1; 58547c6ae99SBarry Smith } 58647c6ae99SBarry Smith 587ce00eea3SSatish Balay ierr = PetscMalloc((Xe-Xs)*(Ye-Ys)*sizeof(PetscInt),&idx);CHKERRQ(ierr); 588ce00eea3SSatish Balay ierr = PetscLogObjectMemory(da,(Xe-Xs)*(Ye-Ys)*sizeof(PetscInt));CHKERRQ(ierr); 58947c6ae99SBarry Smith 590ce00eea3SSatish Balay nn = 0; 59147c6ae99SBarry Smith xbase = bases[rank]; 59247c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 59347c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 594ce00eea3SSatish Balay x_t = lx[n0 % m]; 59547c6ae99SBarry Smith y_t = ly[(n0/m)]; 59647c6ae99SBarry Smith s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x; 59747c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 59847c6ae99SBarry Smith } 599ac119b13SBarry Smith 60047c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 60147c6ae99SBarry Smith x_t = x; 60247c6ae99SBarry Smith y_t = ly[(n1/m)]; 60347c6ae99SBarry Smith s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t; 60447c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 605ac119b13SBarry Smith } else if (by == DMDA_BOUNDARY_MIRROR) { 606ac119b13SBarry Smith for (j=0; j<x; j++) { idx[nn++] = bases[rank] + x*(s_y - i + 1) + j;} 60747c6ae99SBarry Smith } 608ac119b13SBarry Smith 60947c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 610ce00eea3SSatish Balay x_t = lx[n2 % m]; 61147c6ae99SBarry Smith y_t = ly[(n2/m)]; 61247c6ae99SBarry Smith s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t; 61347c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 61447c6ae99SBarry Smith } 61547c6ae99SBarry Smith } 61647c6ae99SBarry Smith 61747c6ae99SBarry Smith for (i=0; i<y; i++) { 61847c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 619ce00eea3SSatish Balay x_t = lx[n3 % m]; 62047c6ae99SBarry Smith /* y_t = y; */ 62147c6ae99SBarry Smith s_t = bases[n3] + (i+1)*x_t - s_x; 62247c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 623ac119b13SBarry Smith } else if (bx == DMDA_BOUNDARY_MIRROR) { 624ac119b13SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = bases[rank] + x*i + s_x - j;} 62547c6ae99SBarry Smith } 62647c6ae99SBarry Smith 62747c6ae99SBarry Smith for (j=0; j<x; j++) { idx[nn++] = xbase++; } /* interior */ 62847c6ae99SBarry Smith 62947c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 630ce00eea3SSatish Balay x_t = lx[n5 % m]; 63147c6ae99SBarry Smith /* y_t = y; */ 63247c6ae99SBarry Smith s_t = bases[n5] + (i)*x_t; 63347c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 634ac119b13SBarry Smith } else if (bx == DMDA_BOUNDARY_MIRROR) { 635ac119b13SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = bases[rank] + x*(i + 1) - 2 - j;} 63647c6ae99SBarry Smith } 63747c6ae99SBarry Smith } 63847c6ae99SBarry Smith 63947c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 64047c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 641ce00eea3SSatish Balay x_t = lx[n6 % m]; 64247c6ae99SBarry Smith /* y_t = ly[(n6/m)]; */ 64347c6ae99SBarry Smith s_t = bases[n6] + (i)*x_t - s_x; 64447c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 64547c6ae99SBarry Smith } 646ac119b13SBarry Smith 64747c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 64847c6ae99SBarry Smith x_t = x; 64947c6ae99SBarry Smith /* y_t = ly[(n7/m)]; */ 65047c6ae99SBarry Smith s_t = bases[n7] + (i-1)*x_t; 65147c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 652ac119b13SBarry Smith } else if (by == DMDA_BOUNDARY_MIRROR){ 653ac119b13SBarry Smith for (j=0; j<x; j++) { idx[nn++] = bases[rank] + x*(y - i - 1) + j;} 65447c6ae99SBarry Smith } 655ac119b13SBarry Smith 65647c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 657ce00eea3SSatish Balay x_t = lx[n8 % m]; 65847c6ae99SBarry Smith /* y_t = ly[(n8/m)]; */ 65947c6ae99SBarry Smith s_t = bases[n8] + (i-1)*x_t; 66047c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 66147c6ae99SBarry Smith } 66247c6ae99SBarry Smith } 66347c6ae99SBarry Smith 664ce00eea3SSatish Balay ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_COPY_VALUES,&from);CHKERRQ(ierr); 66547c6ae99SBarry Smith ierr = VecScatterCreate(global,from,local,to,>ol);CHKERRQ(ierr); 66647c6ae99SBarry Smith ierr = PetscLogObjectParent(da,gtol);CHKERRQ(ierr); 667fcfd50ebSBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 668fcfd50ebSBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 66947c6ae99SBarry Smith 67088661749SPeter Brune if (stencil_type == DMDA_STENCIL_STAR && o == 0) { 671ce00eea3SSatish Balay n0 = sn0; n2 = sn2; n6 = sn6; n8 = sn8; 672ce00eea3SSatish Balay } 673ce00eea3SSatish Balay 67488661749SPeter Brune if (((stencil_type == DMDA_STENCIL_STAR) || 6751321219cSEthan Coon (bx && bx != DMDA_BOUNDARY_PERIODIC) || 67688661749SPeter Brune (by && by != DMDA_BOUNDARY_PERIODIC)) && o == 0) { 67747c6ae99SBarry Smith /* 67847c6ae99SBarry Smith Recompute the local to global mappings, this time keeping the 679ce00eea3SSatish Balay information about the cross corner processor numbers and any ghosted 680ce00eea3SSatish Balay but not periodic indices. 68147c6ae99SBarry Smith */ 68247c6ae99SBarry Smith nn = 0; 68347c6ae99SBarry Smith xbase = bases[rank]; 68447c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 68547c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 686ce00eea3SSatish Balay x_t = lx[n0 % m]; 68747c6ae99SBarry Smith y_t = ly[(n0/m)]; 68847c6ae99SBarry Smith s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x; 68947c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 690ce00eea3SSatish Balay } else if (xs-Xs > 0 && ys-Ys > 0) { 691ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 69247c6ae99SBarry Smith } 69347c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 69447c6ae99SBarry Smith x_t = x; 69547c6ae99SBarry Smith y_t = ly[(n1/m)]; 69647c6ae99SBarry Smith s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t; 69747c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 698ce00eea3SSatish Balay } else if (ys-Ys > 0) { 699624904c4SBarry Smith if (by == DMDA_BOUNDARY_MIRROR) { 700624904c4SBarry Smith for (j=0; j<x; j++) { idx[nn++] = bases[rank] + x*(s_y - i + 1) + j;} 701624904c4SBarry Smith } else { 702ce00eea3SSatish Balay for (j=0; j<x; j++) { idx[nn++] = -1;} 70347c6ae99SBarry Smith } 704624904c4SBarry Smith } 70547c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 706ce00eea3SSatish Balay x_t = lx[n2 % m]; 70747c6ae99SBarry Smith y_t = ly[(n2/m)]; 70847c6ae99SBarry Smith s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t; 70947c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 710ce00eea3SSatish Balay } else if (Xe-xe> 0 && ys-Ys > 0) { 711ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 71247c6ae99SBarry Smith } 71347c6ae99SBarry Smith } 71447c6ae99SBarry Smith 71547c6ae99SBarry Smith for (i=0; i<y; i++) { 71647c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 717ce00eea3SSatish Balay x_t = lx[n3 % m]; 71847c6ae99SBarry Smith /* y_t = y; */ 71947c6ae99SBarry Smith s_t = bases[n3] + (i+1)*x_t - s_x; 72047c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 721ce00eea3SSatish Balay } else if (xs-Xs > 0) { 722624904c4SBarry Smith if (bx == DMDA_BOUNDARY_MIRROR) { 723624904c4SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = bases[rank] + x*i + s_x - j;} 724624904c4SBarry Smith } else { 725ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 72647c6ae99SBarry Smith } 727624904c4SBarry Smith } 72847c6ae99SBarry Smith 72947c6ae99SBarry Smith for (j=0; j<x; j++) { idx[nn++] = xbase++; } /* interior */ 73047c6ae99SBarry Smith 73147c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 732ce00eea3SSatish Balay x_t = lx[n5 % m]; 73347c6ae99SBarry Smith /* y_t = y; */ 73447c6ae99SBarry Smith s_t = bases[n5] + (i)*x_t; 73547c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 736ce00eea3SSatish Balay } else if (Xe-xe > 0) { 737624904c4SBarry Smith if (bx == DMDA_BOUNDARY_MIRROR) { 738624904c4SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = bases[rank] + x*(i + 1) - 2 - j;} 739624904c4SBarry Smith } else { 740ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 74147c6ae99SBarry Smith } 74247c6ae99SBarry Smith } 743624904c4SBarry Smith } 74447c6ae99SBarry Smith 74547c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 74647c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 747ce00eea3SSatish Balay x_t = lx[n6 % m]; 74847c6ae99SBarry Smith /* y_t = ly[(n6/m)]; */ 74947c6ae99SBarry Smith s_t = bases[n6] + (i)*x_t - s_x; 75047c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 751ce00eea3SSatish Balay } else if (xs-Xs > 0 && Ye-ye > 0) { 752ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 75347c6ae99SBarry Smith } 75447c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 75547c6ae99SBarry Smith x_t = x; 75647c6ae99SBarry Smith /* y_t = ly[(n7/m)]; */ 75747c6ae99SBarry Smith s_t = bases[n7] + (i-1)*x_t; 75847c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 759ce00eea3SSatish Balay } else if (Ye-ye > 0) { 760624904c4SBarry Smith if (by == DMDA_BOUNDARY_MIRROR) { 761624904c4SBarry Smith for (j=0; j<x; j++) { idx[nn++] = bases[rank] + x*(y - i - 1) + j;} 762624904c4SBarry Smith } else { 763ce00eea3SSatish Balay for (j=0; j<x; j++) { idx[nn++] = -1;} 76447c6ae99SBarry Smith } 765624904c4SBarry Smith } 76647c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 767ce00eea3SSatish Balay x_t = lx[n8 % m]; 76847c6ae99SBarry Smith /* y_t = ly[(n8/m)]; */ 76947c6ae99SBarry Smith s_t = bases[n8] + (i-1)*x_t; 77047c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 771ce00eea3SSatish Balay } else if (Xe-xe > 0 && Ye-ye > 0) { 772ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 77347c6ae99SBarry Smith } 77447c6ae99SBarry Smith } 77547c6ae99SBarry Smith } 776ce00eea3SSatish Balay /* 777ce00eea3SSatish Balay Set the local to global ordering in the global vector, this allows use 778ce00eea3SSatish Balay of VecSetValuesLocal(). 779ce00eea3SSatish Balay */ 780ce00eea3SSatish Balay ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_OWN_POINTER,<ogis);CHKERRQ(ierr); 781ce00eea3SSatish Balay ierr = PetscMalloc(nn*dof*sizeof(PetscInt),&idx_cpy);CHKERRQ(ierr); 782db87c5ecSEthan Coon ierr = PetscLogObjectMemory(da,nn*dof*sizeof(PetscInt));CHKERRQ(ierr); 783*3bf036e2SBarry Smith ierr = ISGetIndices(ltogis, &idx_full);CHKERRQ(ierr); 784ce00eea3SSatish Balay ierr = PetscMemcpy(idx_cpy,idx_full,nn*dof*sizeof(PetscInt));CHKERRQ(ierr); 785*3bf036e2SBarry Smith ierr = ISRestoreIndices(ltogis, &idx_full);CHKERRQ(ierr); 786ce00eea3SSatish Balay ierr = ISLocalToGlobalMappingCreateIS(ltogis,&da->ltogmap);CHKERRQ(ierr); 787ce00eea3SSatish Balay ierr = PetscLogObjectParent(da,da->ltogmap);CHKERRQ(ierr); 788fcfd50ebSBarry Smith ierr = ISDestroy(<ogis);CHKERRQ(ierr); 789ce00eea3SSatish Balay ierr = ISLocalToGlobalMappingBlock(da->ltogmap,dd->w,&da->ltogmapb);CHKERRQ(ierr); 790ce00eea3SSatish Balay ierr = PetscLogObjectParent(da,da->ltogmap);CHKERRQ(ierr); 79147c6ae99SBarry Smith 792ce00eea3SSatish Balay ierr = PetscFree2(bases,ldims);CHKERRQ(ierr); 79347c6ae99SBarry Smith dd->m = m; dd->n = n; 794ce00eea3SSatish Balay /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */ 795ce00eea3SSatish Balay dd->xs = xs*dof; dd->xe = xe*dof; dd->ys = ys; dd->ye = ye; dd->zs = 0; dd->ze = 1; 796ce00eea3SSatish Balay dd->Xs = Xs*dof; dd->Xe = Xe*dof; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = 0; dd->Ze = 1; 79747c6ae99SBarry Smith 798fcfd50ebSBarry Smith ierr = VecDestroy(&local);CHKERRQ(ierr); 799fcfd50ebSBarry Smith ierr = VecDestroy(&global);CHKERRQ(ierr); 80047c6ae99SBarry Smith 80147c6ae99SBarry Smith dd->gtol = gtol; 80247c6ae99SBarry Smith dd->ltog = ltog; 803ce00eea3SSatish Balay dd->idx = idx_cpy; 804ce00eea3SSatish Balay dd->Nl = nn*dof; 80547c6ae99SBarry Smith dd->base = base; 8069a42bb27SBarry Smith da->ops->view = DMView_DA_2d; 80747c6ae99SBarry Smith dd->ltol = PETSC_NULL; 80847c6ae99SBarry Smith dd->ao = PETSC_NULL; 80947c6ae99SBarry Smith 81047c6ae99SBarry Smith PetscFunctionReturn(0); 81147c6ae99SBarry Smith } 81247c6ae99SBarry Smith 81347c6ae99SBarry Smith #undef __FUNCT__ 814aa219208SBarry Smith #define __FUNCT__ "DMDACreate2d" 81547c6ae99SBarry Smith /*@C 816aa219208SBarry Smith DMDACreate2d - Creates an object that will manage the communication of two-dimensional 81747c6ae99SBarry Smith regular array data that is distributed across some processors. 81847c6ae99SBarry Smith 81947c6ae99SBarry Smith Collective on MPI_Comm 82047c6ae99SBarry Smith 82147c6ae99SBarry Smith Input Parameters: 82247c6ae99SBarry Smith + comm - MPI communicator 8231321219cSEthan Coon . bx,by - type of ghost nodes the array have. 8241321219cSEthan Coon Use one of DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_PERIODIC. 825aa219208SBarry Smith . stencil_type - stencil type. Use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR. 82647c6ae99SBarry 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 82747c6ae99SBarry Smith from the command line with -da_grid_x <M> -da_grid_y <N>) 82847c6ae99SBarry Smith . m,n - corresponding number of processors in each dimension 82947c6ae99SBarry Smith (or PETSC_DECIDE to have calculated) 83047c6ae99SBarry Smith . dof - number of degrees of freedom per node 83147c6ae99SBarry Smith . s - stencil width 83247c6ae99SBarry Smith - lx, ly - arrays containing the number of nodes in each cell along 83347c6ae99SBarry Smith the x and y coordinates, or PETSC_NULL. If non-null, these 83447c6ae99SBarry Smith must be of length as m and n, and the corresponding 83547c6ae99SBarry Smith m and n cannot be PETSC_DECIDE. The sum of the lx[] entries 83647c6ae99SBarry Smith must be M, and the sum of the ly[] entries must be N. 83747c6ae99SBarry Smith 83847c6ae99SBarry Smith Output Parameter: 83947c6ae99SBarry Smith . da - the resulting distributed array object 84047c6ae99SBarry Smith 84147c6ae99SBarry Smith Options Database Key: 842706a11cbSBarry Smith + -dm_view - Calls DMView() at the conclusion of DMDACreate2d() 84347c6ae99SBarry Smith . -da_grid_x <nx> - number of grid points in x direction, if M < 0 84447c6ae99SBarry Smith . -da_grid_y <ny> - number of grid points in y direction, if N < 0 84547c6ae99SBarry Smith . -da_processors_x <nx> - number of processors in x direction 84647c6ae99SBarry Smith . -da_processors_y <ny> - number of processors in y direction 847e0f5d30fSBarry Smith . -da_refine_x <rx> - refinement ratio in x direction 848e0f5d30fSBarry Smith . -da_refine_y <ry> - refinement ratio in y direction 849e0f5d30fSBarry Smith - -da_refine <n> - refine the DMDA n times before creating, if M or N < 0 850e0f5d30fSBarry Smith 85147c6ae99SBarry Smith 85247c6ae99SBarry Smith Level: beginner 85347c6ae99SBarry Smith 85447c6ae99SBarry Smith Notes: 855aa219208SBarry Smith The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the 856aa219208SBarry Smith standard 5-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes 85747c6ae99SBarry Smith the standard 9-pt stencil. 85847c6ae99SBarry Smith 859aa219208SBarry Smith The array data itself is NOT stored in the DMDA, it is stored in Vec objects; 860564755cdSBarry Smith The appropriate vector objects can be obtained with calls to DMCreateGlobalVector() 861564755cdSBarry Smith and DMCreateLocalVector() and calls to VecDuplicate() if more are needed. 86247c6ae99SBarry Smith 86347c6ae99SBarry Smith .keywords: distributed array, create, two-dimensional 86447c6ae99SBarry Smith 865aa219208SBarry Smith .seealso: DMDestroy(), DMView(), DMDACreate1d(), DMDACreate3d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(), 866aa219208SBarry Smith DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMDALocalToLocalBegin(), DMDALocalToLocalEnd(), DMDASetRefinementFactor(), 867d461ba97SBarry Smith DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges() 86847c6ae99SBarry Smith 86947c6ae99SBarry Smith @*/ 870fe16a2e9SBarry Smith 8711321219cSEthan Coon PetscErrorCode DMDACreate2d(MPI_Comm comm,DMDABoundaryType bx,DMDABoundaryType by,DMDAStencilType stencil_type, 8729a42bb27SBarry Smith PetscInt M,PetscInt N,PetscInt m,PetscInt n,PetscInt dof,PetscInt s,const PetscInt lx[],const PetscInt ly[],DM *da) 87347c6ae99SBarry Smith { 87447c6ae99SBarry Smith PetscErrorCode ierr; 87547c6ae99SBarry Smith 87647c6ae99SBarry Smith PetscFunctionBegin; 877aa219208SBarry Smith ierr = DMDACreate(comm, da);CHKERRQ(ierr); 878aa219208SBarry Smith ierr = DMDASetDim(*da, 2);CHKERRQ(ierr); 879aa219208SBarry Smith ierr = DMDASetSizes(*da, M, N, 1);CHKERRQ(ierr); 880aa219208SBarry Smith ierr = DMDASetNumProcs(*da, m, n, PETSC_DECIDE);CHKERRQ(ierr); 881755f258dSLisandro Dalcin ierr = DMDASetBoundaryType(*da, bx, by, DMDA_BOUNDARY_NONE);CHKERRQ(ierr); 882aa219208SBarry Smith ierr = DMDASetDof(*da, dof);CHKERRQ(ierr); 883aa219208SBarry Smith ierr = DMDASetStencilType(*da, stencil_type);CHKERRQ(ierr); 884aa219208SBarry Smith ierr = DMDASetStencilWidth(*da, s);CHKERRQ(ierr); 885aa219208SBarry Smith ierr = DMDASetOwnershipRanges(*da, lx, ly, PETSC_NULL);CHKERRQ(ierr); 88647c6ae99SBarry Smith /* This violates the behavior for other classes, but right now users expect negative dimensions to be handled this way */ 8879a42bb27SBarry Smith ierr = DMSetFromOptions(*da);CHKERRQ(ierr); 8889a42bb27SBarry Smith ierr = DMSetUp(*da);CHKERRQ(ierr); 889ca266f36SBarry Smith ierr = DMViewFromOptions(*da,"-dm_view");CHKERRQ(ierr); 89047c6ae99SBarry Smith PetscFunctionReturn(0); 89147c6ae99SBarry Smith } 892