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 PetscMPIInt rank; 8c9493c35SStefano Zampini PetscBool iascii,isdraw,isglvis,isbinary; 947c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 109a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 119a42bb27SBarry Smith PetscBool ismatlab; 129a42bb27SBarry Smith #endif 1347c6ae99SBarry Smith 1447c6ae99SBarry Smith PetscFunctionBegin; 159566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)da),&rank)); 1647c6ae99SBarry Smith 179566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 189566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw)); 199566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERGLVIS,&isglvis)); 209566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary)); 219a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 229566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERMATLAB,&ismatlab)); 239a42bb27SBarry Smith #endif 2447c6ae99SBarry Smith if (iascii) { 2547c6ae99SBarry Smith PetscViewerFormat format; 2647c6ae99SBarry Smith 279566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format)); 2876a8abe0SBarry Smith if (format == PETSC_VIEWER_LOAD_BALANCE) { 2976a8abe0SBarry Smith PetscInt i,nmax = 0,nmin = PETSC_MAX_INT,navg = 0,*nz,nzlocal; 3076a8abe0SBarry Smith DMDALocalInfo info; 3176a8abe0SBarry Smith PetscMPIInt size; 329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)da),&size)); 339566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da,&info)); 3476a8abe0SBarry Smith nzlocal = info.xm*info.ym; 359566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size,&nz)); 369566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&nzlocal,1,MPIU_INT,nz,1,MPIU_INT,PetscObjectComm((PetscObject)da))); 3776a8abe0SBarry Smith for (i=0; i<(PetscInt)size; i++) { 3876a8abe0SBarry Smith nmax = PetscMax(nmax,nz[i]); 3976a8abe0SBarry Smith nmin = PetscMin(nmin,nz[i]); 4076a8abe0SBarry Smith navg += nz[i]; 4176a8abe0SBarry Smith } 429566063dSJacob Faibussowitsch PetscCall(PetscFree(nz)); 4376a8abe0SBarry Smith navg = navg/size; 4463a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Load Balance - Grid Points: Min %" PetscInt_FMT " avg %" PetscInt_FMT " max %" PetscInt_FMT "\n",nmin,navg,nmax)); 4576a8abe0SBarry Smith PetscFunctionReturn(0); 4676a8abe0SBarry Smith } 478ec8862eSJed Brown if (format != PETSC_VIEWER_ASCII_VTK_DEPRECATED && format != PETSC_VIEWER_ASCII_VTK_CELL_DEPRECATED && format != PETSC_VIEWER_ASCII_GLVIS) { 48aa219208SBarry Smith DMDALocalInfo info; 499566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da,&info)); 509566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 5163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"Processor [%d] M %" PetscInt_FMT " N %" PetscInt_FMT " m %" PetscInt_FMT " n %" PetscInt_FMT " w %" PetscInt_FMT " s %" PetscInt_FMT "\n",rank,dd->M,dd->N,dd->m,dd->n,dd->w,dd->s)); 5263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"X range of indices: %" PetscInt_FMT " %" PetscInt_FMT ", Y range of indices: %" PetscInt_FMT " %" PetscInt_FMT "\n",info.xs,info.xs+info.xm,info.ys,info.ys+info.ym)); 539566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 549566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 55*1baa6e33SBarry Smith } else if (format == PETSC_VIEWER_ASCII_GLVIS) PetscCall(DMView_DA_GLVis(da,viewer)); 56*1baa6e33SBarry Smith else PetscCall(DMView_DA_VTK(da,viewer)); 5747c6ae99SBarry Smith } else if (isdraw) { 5847c6ae99SBarry Smith PetscDraw draw; 5947c6ae99SBarry Smith double ymin = -1*dd->s-1,ymax = dd->N+dd->s; 6047c6ae99SBarry Smith double xmin = -1*dd->s-1,xmax = dd->M+dd->s; 6147c6ae99SBarry Smith double x,y; 628ea3bf28SBarry Smith PetscInt base; 638ea3bf28SBarry Smith const PetscInt *idx; 6447c6ae99SBarry Smith char node[10]; 6547c6ae99SBarry Smith PetscBool isnull; 6647c6ae99SBarry Smith 679566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(viewer,0,&draw)); 689566063dSJacob Faibussowitsch PetscCall(PetscDrawIsNull(draw,&isnull)); 695b399a63SLisandro Dalcin if (isnull) PetscFunctionReturn(0); 7047c6ae99SBarry Smith 719566063dSJacob Faibussowitsch PetscCall(PetscDrawCheckResizedWindow(draw)); 729566063dSJacob Faibussowitsch PetscCall(PetscDrawClear(draw)); 739566063dSJacob Faibussowitsch PetscCall(PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax)); 745b399a63SLisandro Dalcin 75d0609cedSBarry Smith PetscDrawCollectiveBegin(draw); 7647c6ae99SBarry Smith /* first processor draw all node lines */ 77dd400576SPatrick Sanan if (rank == 0) { 7847c6ae99SBarry Smith ymin = 0.0; ymax = dd->N - 1; 7947c6ae99SBarry Smith for (xmin=0; xmin<dd->M; xmin++) { 809566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK)); 8147c6ae99SBarry Smith } 8247c6ae99SBarry Smith xmin = 0.0; xmax = dd->M - 1; 8347c6ae99SBarry Smith for (ymin=0; ymin<dd->N; ymin++) { 849566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK)); 8547c6ae99SBarry Smith } 8647c6ae99SBarry Smith } 87d0609cedSBarry Smith PetscDrawCollectiveEnd(draw); 889566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 899566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw)); 9047c6ae99SBarry Smith 91d0609cedSBarry Smith PetscDrawCollectiveBegin(draw); 9247c6ae99SBarry Smith /* draw my box */ 935b399a63SLisandro Dalcin xmin = dd->xs/dd->w; xmax =(dd->xe-1)/dd->w; ymin = dd->ys; ymax = dd->ye - 1; 949566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED)); 959566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED)); 969566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED)); 979566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED)); 9847c6ae99SBarry Smith /* put in numbers */ 9947c6ae99SBarry Smith base = (dd->base)/dd->w; 10047c6ae99SBarry Smith for (y=ymin; y<=ymax; y++) { 10147c6ae99SBarry Smith for (x=xmin; x<=xmax; x++) { 1029566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(node,sizeof(node),"%d",(int)base++)); 1039566063dSJacob Faibussowitsch PetscCall(PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node)); 10447c6ae99SBarry Smith } 10547c6ae99SBarry Smith } 106d0609cedSBarry Smith PetscDrawCollectiveEnd(draw); 1079566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 1089566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw)); 10947c6ae99SBarry Smith 110d0609cedSBarry Smith PetscDrawCollectiveBegin(draw); 1115b399a63SLisandro Dalcin /* overlay ghost numbers, useful for error checking */ 1129566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(da->ltogmap,&idx)); 1135b399a63SLisandro Dalcin base = 0; xmin = dd->Xs; xmax = dd->Xe; ymin = dd->Ys; ymax = dd->Ye; 11447c6ae99SBarry Smith for (y=ymin; y<ymax; y++) { 11547c6ae99SBarry Smith for (x=xmin; x<xmax; x++) { 11647c6ae99SBarry Smith if ((base % dd->w) == 0) { 1179566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(node,sizeof(node),"%d",(int)(idx[base/dd->w]))); 1189566063dSJacob Faibussowitsch PetscCall(PetscDrawString(draw,x/dd->w,y,PETSC_DRAW_BLUE,node)); 11947c6ae99SBarry Smith } 12047c6ae99SBarry Smith base++; 12147c6ae99SBarry Smith } 12247c6ae99SBarry Smith } 1239566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(da->ltogmap,&idx)); 124d0609cedSBarry Smith PetscDrawCollectiveEnd(draw); 1259566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 1269566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw)); 1279566063dSJacob Faibussowitsch PetscCall(PetscDrawSave(draw)); 128c9493c35SStefano Zampini } else if (isglvis) { 1299566063dSJacob Faibussowitsch PetscCall(DMView_DA_GLVis(da,viewer)); 1309a42bb27SBarry Smith } else if (isbinary) { 1319566063dSJacob Faibussowitsch PetscCall(DMView_DA_Binary(da,viewer)); 1329a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 1339a42bb27SBarry Smith } else if (ismatlab) { 1349566063dSJacob Faibussowitsch PetscCall(DMView_DA_Matlab(da,viewer)); 1359a42bb27SBarry Smith #endif 13611aeaf0aSBarry Smith } 13747c6ae99SBarry Smith PetscFunctionReturn(0); 13847c6ae99SBarry Smith } 13947c6ae99SBarry Smith 14047c6ae99SBarry Smith #if defined(new) 14147c6ae99SBarry Smith /* 142aa219208SBarry Smith DMDAGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix where local 143aa219208SBarry Smith function lives on a DMDA 14447c6ae99SBarry Smith 14547c6ae99SBarry Smith y ~= (F(u + ha) - F(u))/h, 14647c6ae99SBarry Smith where F = nonlinear function, as set by SNESSetFunction() 14747c6ae99SBarry Smith u = current iterate 14847c6ae99SBarry Smith h = difference interval 14947c6ae99SBarry Smith */ 150aa219208SBarry Smith PetscErrorCode DMDAGetDiagonal_MFFD(DM da,Vec U,Vec a) 15147c6ae99SBarry Smith { 15247c6ae99SBarry Smith PetscScalar h,*aa,*ww,v; 15347c6ae99SBarry Smith PetscReal epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON; 15447c6ae99SBarry Smith PetscInt gI,nI; 15547c6ae99SBarry Smith MatStencil stencil; 156aa219208SBarry Smith DMDALocalInfo info; 15747c6ae99SBarry Smith 15847c6ae99SBarry Smith PetscFunctionBegin; 1599566063dSJacob Faibussowitsch PetscCall((*ctx->func)(0,U,a,ctx->funcctx)); 1609566063dSJacob Faibussowitsch PetscCall((*ctx->funcisetbase)(U,ctx->funcctx)); 16147c6ae99SBarry Smith 1629566063dSJacob Faibussowitsch PetscCall(VecGetArray(U,&ww)); 1639566063dSJacob Faibussowitsch PetscCall(VecGetArray(a,&aa)); 16447c6ae99SBarry Smith 16547c6ae99SBarry Smith nI = 0; 16647c6ae99SBarry Smith h = ww[gI]; 16747c6ae99SBarry Smith if (h == 0.0) h = 1.0; 16847c6ae99SBarry Smith if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin; 16947c6ae99SBarry Smith else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin; 17047c6ae99SBarry Smith h *= epsilon; 17147c6ae99SBarry Smith 17247c6ae99SBarry Smith ww[gI] += h; 1739566063dSJacob Faibussowitsch PetscCall((*ctx->funci)(i,w,&v,ctx->funcctx)); 17447c6ae99SBarry Smith aa[nI] = (v - aa[nI])/h; 17547c6ae99SBarry Smith ww[gI] -= h; 17647c6ae99SBarry Smith nI++; 1778865f1eaSKarl Rupp 1789566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U,&ww)); 1799566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(a,&aa)); 18047c6ae99SBarry Smith PetscFunctionReturn(0); 18147c6ae99SBarry Smith } 18247c6ae99SBarry Smith #endif 18347c6ae99SBarry Smith 1847087cfbeSBarry Smith PetscErrorCode DMSetUp_DA_2D(DM da) 18547c6ae99SBarry Smith { 18647c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 18747c6ae99SBarry Smith const PetscInt M = dd->M; 18847c6ae99SBarry Smith const PetscInt N = dd->N; 18947c6ae99SBarry Smith PetscInt m = dd->m; 19047c6ae99SBarry Smith PetscInt n = dd->n; 19147c6ae99SBarry Smith const PetscInt dof = dd->w; 19247c6ae99SBarry Smith const PetscInt s = dd->s; 193bff4a2f0SMatthew G. Knepley DMBoundaryType bx = dd->bx; 194bff4a2f0SMatthew G. Knepley DMBoundaryType by = dd->by; 19519fd82e9SBarry Smith DMDAStencilType stencil_type = dd->stencil_type; 19647c6ae99SBarry Smith PetscInt *lx = dd->lx; 19747c6ae99SBarry Smith PetscInt *ly = dd->ly; 19847c6ae99SBarry Smith MPI_Comm comm; 19947c6ae99SBarry Smith PetscMPIInt rank,size; 200bd1fc5aeSBarry Smith PetscInt xs,xe,ys,ye,x,y,Xs,Xe,Ys,Ye,IXs,IXe,IYs,IYe; 2018ea3bf28SBarry Smith PetscInt up,down,left,right,i,n0,n1,n2,n3,n5,n6,n7,n8,*idx,nn; 202db87c5ecSEthan Coon PetscInt xbase,*bases,*ldims,j,x_t,y_t,s_t,base,count; 20347c6ae99SBarry Smith PetscInt s_x,s_y; /* s proportionalized to w */ 20447c6ae99SBarry Smith PetscInt sn0 = 0,sn2 = 0,sn6 = 0,sn8 = 0; 20547c6ae99SBarry Smith Vec local,global; 206bd1fc5aeSBarry Smith VecScatter gtol; 20745b6f7e9SBarry Smith IS to,from; 20847c6ae99SBarry Smith 20947c6ae99SBarry Smith PetscFunctionBegin; 2107a8be351SBarry Smith PetscCheck(stencil_type != DMDA_STENCIL_BOX || (bx != DM_BOUNDARY_MIRROR && by != DM_BOUNDARY_MIRROR),PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Mirror boundary and box stencil"); 2119566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)da,&comm)); 2123855c12bSBarry Smith #if !defined(PETSC_USE_64BIT_INDICES) 21363a3b9bcSJacob Faibussowitsch PetscCheck(((PetscInt64) M)*((PetscInt64) N)*((PetscInt64) dof) <= (PetscInt64) PETSC_MPI_INT_MAX,comm,PETSC_ERR_INT_OVERFLOW,"Mesh of %" PetscInt_FMT " by %" PetscInt_FMT " by %" PetscInt_FMT " (dof) is too large for 32 bit indices",M,N,dof); 2143855c12bSBarry Smith #endif 21547c6ae99SBarry Smith 2169566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 2179566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm,&rank)); 21847c6ae99SBarry Smith 2197d310018SBarry Smith dd->p = 1; 22047c6ae99SBarry Smith if (m != PETSC_DECIDE) { 22163a3b9bcSJacob Faibussowitsch PetscCheck(m >= 1,comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %" PetscInt_FMT,m); 22263a3b9bcSJacob Faibussowitsch else PetscCheck(m <= size,comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %" PetscInt_FMT " %d",m,size); 22347c6ae99SBarry Smith } 22447c6ae99SBarry Smith if (n != PETSC_DECIDE) { 22563a3b9bcSJacob Faibussowitsch PetscCheck(n >= 1,comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %" PetscInt_FMT,n); 22663a3b9bcSJacob Faibussowitsch else PetscCheck(n <= size,comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %" PetscInt_FMT " %d",n,size); 22747c6ae99SBarry Smith } 22847c6ae99SBarry Smith 22947c6ae99SBarry Smith if (m == PETSC_DECIDE || n == PETSC_DECIDE) { 23047c6ae99SBarry Smith if (n != PETSC_DECIDE) { 23147c6ae99SBarry Smith m = size/n; 23247c6ae99SBarry Smith } else if (m != PETSC_DECIDE) { 23347c6ae99SBarry Smith n = size/m; 23447c6ae99SBarry Smith } else { 23547c6ae99SBarry Smith /* try for squarish distribution */ 236369cc0aeSBarry Smith m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)N))); 23747c6ae99SBarry Smith if (!m) m = 1; 23847c6ae99SBarry Smith while (m > 0) { 23947c6ae99SBarry Smith n = size/m; 24047c6ae99SBarry Smith if (m*n == size) break; 24147c6ae99SBarry Smith m--; 24247c6ae99SBarry Smith } 24347c6ae99SBarry Smith if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;} 24447c6ae99SBarry Smith } 2457a8be351SBarry Smith PetscCheck(m*n == size,comm,PETSC_ERR_PLIB,"Unable to create partition, check the size of the communicator and input m and n "); 2467a8be351SBarry Smith } else PetscCheck(m*n == size,comm,PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition"); 24747c6ae99SBarry Smith 24863a3b9bcSJacob Faibussowitsch PetscCheck(M >= m,comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %" PetscInt_FMT " %" PetscInt_FMT,M,m); 24963a3b9bcSJacob Faibussowitsch PetscCheck(N >= n,comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %" PetscInt_FMT " %" PetscInt_FMT,N,n); 25047c6ae99SBarry Smith 25147c6ae99SBarry Smith /* 25247c6ae99SBarry Smith Determine locally owned region 25347c6ae99SBarry Smith xs is the first local node number, x is the number of local nodes 25447c6ae99SBarry Smith */ 25547c6ae99SBarry Smith if (!lx) { 2569566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m, &dd->lx)); 25747c6ae99SBarry Smith lx = dd->lx; 25847c6ae99SBarry Smith for (i=0; i<m; i++) { 25947c6ae99SBarry Smith lx[i] = M/m + ((M % m) > i); 26047c6ae99SBarry Smith } 26147c6ae99SBarry Smith } 26247c6ae99SBarry Smith x = lx[rank % m]; 26347c6ae99SBarry Smith xs = 0; 26447c6ae99SBarry Smith for (i=0; i<(rank % m); i++) { 26547c6ae99SBarry Smith xs += lx[i]; 26647c6ae99SBarry Smith } 26776bd3646SJed Brown if (PetscDefined(USE_DEBUG)) { 26847c6ae99SBarry Smith left = xs; 26947c6ae99SBarry Smith for (i=(rank % m); i<m; i++) { 27047c6ae99SBarry Smith left += lx[i]; 27147c6ae99SBarry Smith } 27263a3b9bcSJacob Faibussowitsch PetscCheck(left == M,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Sum of lx across processors not equal to M: %" PetscInt_FMT " %" PetscInt_FMT,left,M); 27376bd3646SJed Brown } 27447c6ae99SBarry Smith 27547c6ae99SBarry Smith /* 27647c6ae99SBarry Smith Determine locally owned region 27747c6ae99SBarry Smith ys is the first local node number, y is the number of local nodes 27847c6ae99SBarry Smith */ 27947c6ae99SBarry Smith if (!ly) { 2809566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &dd->ly)); 28147c6ae99SBarry Smith ly = dd->ly; 28247c6ae99SBarry Smith for (i=0; i<n; i++) { 28347c6ae99SBarry Smith ly[i] = N/n + ((N % n) > i); 28447c6ae99SBarry Smith } 28547c6ae99SBarry Smith } 28647c6ae99SBarry Smith y = ly[rank/m]; 28747c6ae99SBarry Smith ys = 0; 28847c6ae99SBarry Smith for (i=0; i<(rank/m); i++) { 28947c6ae99SBarry Smith ys += ly[i]; 29047c6ae99SBarry Smith } 29176bd3646SJed Brown if (PetscDefined(USE_DEBUG)) { 29247c6ae99SBarry Smith left = ys; 29347c6ae99SBarry Smith for (i=(rank/m); i<n; i++) { 29447c6ae99SBarry Smith left += ly[i]; 29547c6ae99SBarry Smith } 29663a3b9bcSJacob Faibussowitsch PetscCheck(left == N,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Sum of ly across processors not equal to N: %" PetscInt_FMT " %" PetscInt_FMT,left,N); 29776bd3646SJed Brown } 29847c6ae99SBarry Smith 299bcea557cSEthan Coon /* 300bcea557cSEthan Coon check if the scatter requires more than one process neighbor or wraps around 301bcea557cSEthan Coon the domain more than once 302bcea557cSEthan Coon */ 30363a3b9bcSJacob Faibussowitsch PetscCheck((x >= s) || ((m <= 1) && (bx != DM_BOUNDARY_PERIODIC)),PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local x-width of domain x %" PetscInt_FMT " is smaller than stencil width s %" PetscInt_FMT,x,s); 30463a3b9bcSJacob Faibussowitsch PetscCheck((y >= s) || ((n <= 1) && (by != DM_BOUNDARY_PERIODIC)),PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local y-width of domain y %" PetscInt_FMT " is smaller than stencil width s %" PetscInt_FMT,y,s); 30547c6ae99SBarry Smith xe = xs + x; 30647c6ae99SBarry Smith ye = ys + y; 30747c6ae99SBarry Smith 308ce00eea3SSatish Balay /* determine ghost region (Xs) and region scattered into (IXs) */ 309d9c9ebe5SPeter Brune if (xs-s > 0) { 310d9c9ebe5SPeter Brune Xs = xs - s; IXs = xs - s; 31188661749SPeter Brune } else { 31288661749SPeter Brune if (bx) { 31388661749SPeter Brune Xs = xs - s; 31488661749SPeter Brune } else { 31588661749SPeter Brune Xs = 0; 31688661749SPeter Brune } 31788661749SPeter Brune IXs = 0; 31888661749SPeter Brune } 319d9c9ebe5SPeter Brune if (xe+s <= M) { 320d9c9ebe5SPeter Brune Xe = xe + s; IXe = xe + s; 32188661749SPeter Brune } else { 32288661749SPeter Brune if (bx) { 323d9c9ebe5SPeter Brune Xs = xs - s; Xe = xe + s; 32488661749SPeter Brune } else { 32588661749SPeter Brune Xe = M; 32688661749SPeter Brune } 32788661749SPeter Brune IXe = M; 32888661749SPeter Brune } 32947c6ae99SBarry Smith 330bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC || bx == DM_BOUNDARY_MIRROR) { 331d9c9ebe5SPeter Brune IXs = xs - s; 332d9c9ebe5SPeter Brune IXe = xe + s; 333d9c9ebe5SPeter Brune Xs = xs - s; 334d9c9ebe5SPeter Brune Xe = xe + s; 33588661749SPeter Brune } 33647c6ae99SBarry Smith 337d9c9ebe5SPeter Brune if (ys-s > 0) { 338d9c9ebe5SPeter Brune Ys = ys - s; IYs = ys - s; 33988661749SPeter Brune } else { 34088661749SPeter Brune if (by) { 34188661749SPeter Brune Ys = ys - s; 34288661749SPeter Brune } else { 34388661749SPeter Brune Ys = 0; 34488661749SPeter Brune } 34588661749SPeter Brune IYs = 0; 34688661749SPeter Brune } 347d9c9ebe5SPeter Brune if (ye+s <= N) { 348d9c9ebe5SPeter Brune Ye = ye + s; IYe = ye + s; 34988661749SPeter Brune } else { 35088661749SPeter Brune if (by) { 35188661749SPeter Brune Ye = ye + s; 35288661749SPeter Brune } else { 35388661749SPeter Brune Ye = N; 35488661749SPeter Brune } 35588661749SPeter Brune IYe = N; 35688661749SPeter Brune } 35788661749SPeter Brune 358bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_PERIODIC || by == DM_BOUNDARY_MIRROR) { 359d9c9ebe5SPeter Brune IYs = ys - s; 360d9c9ebe5SPeter Brune IYe = ye + s; 361d9c9ebe5SPeter Brune Ys = ys - s; 362d9c9ebe5SPeter Brune Ye = ye + s; 36388661749SPeter Brune } 36488661749SPeter Brune 36588661749SPeter Brune /* stencil length in each direction */ 366d9c9ebe5SPeter Brune s_x = s; 367d9c9ebe5SPeter Brune s_y = s; 36847c6ae99SBarry Smith 36947c6ae99SBarry Smith /* determine starting point of each processor */ 37047c6ae99SBarry Smith nn = x*y; 3719566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(size+1,&bases,size,&ldims)); 3729566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm)); 37347c6ae99SBarry Smith bases[0] = 0; 37447c6ae99SBarry Smith for (i=1; i<=size; i++) { 37547c6ae99SBarry Smith bases[i] = ldims[i-1]; 37647c6ae99SBarry Smith } 37747c6ae99SBarry Smith for (i=1; i<=size; i++) { 37847c6ae99SBarry Smith bases[i] += bases[i-1]; 37947c6ae99SBarry Smith } 380ce00eea3SSatish Balay base = bases[rank]*dof; 38147c6ae99SBarry Smith 38247c6ae99SBarry Smith /* allocate the base parallel and sequential vectors */ 383ce00eea3SSatish Balay dd->Nlocal = x*y*dof; 3849566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,NULL,&global)); 385ce00eea3SSatish Balay dd->nlocal = (Xe-Xs)*(Ye-Ys)*dof; 3869566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,NULL,&local)); 38747c6ae99SBarry Smith 3884104a7a0SPatrick Sanan /* generate global to local vector scatter and local to global mapping*/ 38947c6ae99SBarry Smith 390ce00eea3SSatish Balay /* global to local must include ghost points within the domain, 391ce00eea3SSatish Balay but not ghost points outside the domain that aren't periodic */ 3929566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((IXe-IXs)*(IYe-IYs),&idx)); 393d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_BOX) { 394ce00eea3SSatish Balay left = IXs - Xs; right = left + (IXe-IXs); 395ce00eea3SSatish Balay down = IYs - Ys; up = down + (IYe-IYs); 396ce00eea3SSatish Balay count = 0; 397ce00eea3SSatish Balay for (i=down; i<up; i++) { 398ce00eea3SSatish Balay for (j=left; j<right; j++) { 399ce00eea3SSatish Balay idx[count++] = j + i*(Xe-Xs); 400ce00eea3SSatish Balay } 401ce00eea3SSatish Balay } 4029566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to)); 403ce00eea3SSatish Balay 40447c6ae99SBarry Smith } else { 40547c6ae99SBarry Smith /* must drop into cross shape region */ 40647c6ae99SBarry Smith /* ---------| 40747c6ae99SBarry Smith | top | 408ce00eea3SSatish Balay |--- ---| up 40947c6ae99SBarry Smith | middle | 41047c6ae99SBarry Smith | | 411ce00eea3SSatish Balay ---- ---- down 41247c6ae99SBarry Smith | bottom | 41347c6ae99SBarry Smith ----------- 41447c6ae99SBarry Smith Xs xs xe Xe */ 415ce00eea3SSatish Balay left = xs - Xs; right = left + x; 416ce00eea3SSatish Balay down = ys - Ys; up = down + y; 41747c6ae99SBarry Smith count = 0; 418ce00eea3SSatish Balay /* bottom */ 419ce00eea3SSatish Balay for (i=(IYs-Ys); i<down; i++) { 420ce00eea3SSatish Balay for (j=left; j<right; j++) { 421ce00eea3SSatish Balay idx[count++] = j + i*(Xe-Xs); 42247c6ae99SBarry Smith } 42347c6ae99SBarry Smith } 42447c6ae99SBarry Smith /* middle */ 42547c6ae99SBarry Smith for (i=down; i<up; i++) { 426ce00eea3SSatish Balay for (j=(IXs-Xs); j<(IXe-Xs); j++) { 427ce00eea3SSatish Balay idx[count++] = j + i*(Xe-Xs); 42847c6ae99SBarry Smith } 42947c6ae99SBarry Smith } 43047c6ae99SBarry Smith /* top */ 431ce00eea3SSatish Balay for (i=up; i<up+IYe-ye; i++) { 432ce00eea3SSatish Balay for (j=left; j<right; j++) { 433ce00eea3SSatish Balay idx[count++] = j + i*(Xe-Xs); 43447c6ae99SBarry Smith } 43547c6ae99SBarry Smith } 4369566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to)); 43747c6ae99SBarry Smith } 43847c6ae99SBarry Smith 43947c6ae99SBarry Smith /* determine who lies on each side of us stored in n6 n7 n8 44047c6ae99SBarry Smith n3 n5 44147c6ae99SBarry Smith n0 n1 n2 44247c6ae99SBarry Smith */ 44347c6ae99SBarry Smith 44447c6ae99SBarry Smith /* Assume the Non-Periodic Case */ 44547c6ae99SBarry Smith n1 = rank - m; 44647c6ae99SBarry Smith if (rank % m) { 44747c6ae99SBarry Smith n0 = n1 - 1; 44847c6ae99SBarry Smith } else { 44947c6ae99SBarry Smith n0 = -1; 45047c6ae99SBarry Smith } 45147c6ae99SBarry Smith if ((rank+1) % m) { 45247c6ae99SBarry Smith n2 = n1 + 1; 45347c6ae99SBarry Smith n5 = rank + 1; 45447c6ae99SBarry Smith n8 = rank + m + 1; if (n8 >= m*n) n8 = -1; 45547c6ae99SBarry Smith } else { 45647c6ae99SBarry Smith n2 = -1; n5 = -1; n8 = -1; 45747c6ae99SBarry Smith } 45847c6ae99SBarry Smith if (rank % m) { 45947c6ae99SBarry Smith n3 = rank - 1; 46047c6ae99SBarry Smith n6 = n3 + m; if (n6 >= m*n) n6 = -1; 46147c6ae99SBarry Smith } else { 46247c6ae99SBarry Smith n3 = -1; n6 = -1; 46347c6ae99SBarry Smith } 46447c6ae99SBarry Smith n7 = rank + m; if (n7 >= m*n) n7 = -1; 46547c6ae99SBarry Smith 466bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC && by == DM_BOUNDARY_PERIODIC) { 46747c6ae99SBarry Smith /* Modify for Periodic Cases */ 46847c6ae99SBarry Smith /* Handle all four corners */ 46947c6ae99SBarry Smith if ((n6 < 0) && (n7 < 0) && (n3 < 0)) n6 = m-1; 47047c6ae99SBarry Smith if ((n8 < 0) && (n7 < 0) && (n5 < 0)) n8 = 0; 47147c6ae99SBarry Smith if ((n2 < 0) && (n5 < 0) && (n1 < 0)) n2 = size-m; 47247c6ae99SBarry Smith if ((n0 < 0) && (n3 < 0) && (n1 < 0)) n0 = size-1; 47347c6ae99SBarry Smith 47447c6ae99SBarry Smith /* Handle Top and Bottom Sides */ 47547c6ae99SBarry Smith if (n1 < 0) n1 = rank + m * (n-1); 47647c6ae99SBarry Smith if (n7 < 0) n7 = rank - m * (n-1); 47747c6ae99SBarry Smith if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1; 47847c6ae99SBarry Smith if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1; 47947c6ae99SBarry Smith if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1; 48047c6ae99SBarry Smith if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1; 48147c6ae99SBarry Smith 48247c6ae99SBarry Smith /* Handle Left and Right Sides */ 48347c6ae99SBarry Smith if (n3 < 0) n3 = rank + (m-1); 48447c6ae99SBarry Smith if (n5 < 0) n5 = rank - (m-1); 48547c6ae99SBarry Smith if ((n1 >= 0) && (n0 < 0)) n0 = rank-1; 48647c6ae99SBarry Smith if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1; 48747c6ae99SBarry Smith if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1; 48847c6ae99SBarry Smith if ((n7 >= 0) && (n8 < 0)) n8 = rank+1; 489bff4a2f0SMatthew G. Knepley } else if (by == DM_BOUNDARY_PERIODIC) { /* Handle Top and Bottom Sides */ 490ce00eea3SSatish Balay if (n1 < 0) n1 = rank + m * (n-1); 491ce00eea3SSatish Balay if (n7 < 0) n7 = rank - m * (n-1); 492ce00eea3SSatish Balay if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1; 493ce00eea3SSatish Balay if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1; 494ce00eea3SSatish Balay if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1; 495ce00eea3SSatish Balay if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1; 496bff4a2f0SMatthew G. Knepley } else if (bx == DM_BOUNDARY_PERIODIC) { /* Handle Left and Right Sides */ 497ce00eea3SSatish Balay if (n3 < 0) n3 = rank + (m-1); 498ce00eea3SSatish Balay if (n5 < 0) n5 = rank - (m-1); 499ce00eea3SSatish Balay if ((n1 >= 0) && (n0 < 0)) n0 = rank-1; 500ce00eea3SSatish Balay if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1; 501ce00eea3SSatish Balay if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1; 502ce00eea3SSatish Balay if ((n7 >= 0) && (n8 < 0)) n8 = rank+1; 50347c6ae99SBarry Smith } 504ce00eea3SSatish Balay 5059566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(9,&dd->neighbors)); 5068865f1eaSKarl Rupp 50747c6ae99SBarry Smith dd->neighbors[0] = n0; 50847c6ae99SBarry Smith dd->neighbors[1] = n1; 50947c6ae99SBarry Smith dd->neighbors[2] = n2; 51047c6ae99SBarry Smith dd->neighbors[3] = n3; 51147c6ae99SBarry Smith dd->neighbors[4] = rank; 51247c6ae99SBarry Smith dd->neighbors[5] = n5; 51347c6ae99SBarry Smith dd->neighbors[6] = n6; 51447c6ae99SBarry Smith dd->neighbors[7] = n7; 51547c6ae99SBarry Smith dd->neighbors[8] = n8; 51647c6ae99SBarry Smith 517d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_STAR) { 51847c6ae99SBarry Smith /* save corner processor numbers */ 51947c6ae99SBarry Smith sn0 = n0; sn2 = n2; sn6 = n6; sn8 = n8; 52047c6ae99SBarry Smith n0 = n2 = n6 = n8 = -1; 52147c6ae99SBarry Smith } 52247c6ae99SBarry Smith 5239566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((Xe-Xs)*(Ye-Ys),&idx)); 52447c6ae99SBarry Smith 525ce00eea3SSatish Balay nn = 0; 52647c6ae99SBarry Smith xbase = bases[rank]; 52747c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 52847c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 529ce00eea3SSatish Balay x_t = lx[n0 % m]; 53047c6ae99SBarry Smith y_t = ly[(n0/m)]; 53147c6ae99SBarry Smith s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x; 5328865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 53347c6ae99SBarry Smith } 534ac119b13SBarry Smith 53547c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 53647c6ae99SBarry Smith x_t = x; 53747c6ae99SBarry Smith y_t = ly[(n1/m)]; 53847c6ae99SBarry Smith s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t; 5398865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 540bff4a2f0SMatthew G. Knepley } else if (by == DM_BOUNDARY_MIRROR) { 5418865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(s_y - i + 1) + j; 54247c6ae99SBarry Smith } 543ac119b13SBarry Smith 54447c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 545ce00eea3SSatish Balay x_t = lx[n2 % m]; 54647c6ae99SBarry Smith y_t = ly[(n2/m)]; 54747c6ae99SBarry Smith s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t; 5488865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 54947c6ae99SBarry Smith } 55047c6ae99SBarry Smith } 55147c6ae99SBarry Smith 55247c6ae99SBarry Smith for (i=0; i<y; i++) { 55347c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 554ce00eea3SSatish Balay x_t = lx[n3 % m]; 55547c6ae99SBarry Smith /* y_t = y; */ 55647c6ae99SBarry Smith s_t = bases[n3] + (i+1)*x_t - s_x; 5578865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 558bff4a2f0SMatthew G. Knepley } else if (bx == DM_BOUNDARY_MIRROR) { 5598865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*i + s_x - j; 56047c6ae99SBarry Smith } 56147c6ae99SBarry Smith 5628865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = xbase++; /* interior */ 56347c6ae99SBarry Smith 56447c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 565ce00eea3SSatish Balay x_t = lx[n5 % m]; 56647c6ae99SBarry Smith /* y_t = y; */ 56747c6ae99SBarry Smith s_t = bases[n5] + (i)*x_t; 5688865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 569bff4a2f0SMatthew G. Knepley } else if (bx == DM_BOUNDARY_MIRROR) { 5708865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*(i + 1) - 2 - j; 57147c6ae99SBarry Smith } 57247c6ae99SBarry Smith } 57347c6ae99SBarry Smith 57447c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 57547c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 576ce00eea3SSatish Balay x_t = lx[n6 % m]; 57747c6ae99SBarry Smith /* y_t = ly[(n6/m)]; */ 57847c6ae99SBarry Smith s_t = bases[n6] + (i)*x_t - s_x; 5798865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 58047c6ae99SBarry Smith } 581ac119b13SBarry Smith 58247c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 58347c6ae99SBarry Smith x_t = x; 58447c6ae99SBarry Smith /* y_t = ly[(n7/m)]; */ 58547c6ae99SBarry Smith s_t = bases[n7] + (i-1)*x_t; 5868865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 587bff4a2f0SMatthew G. Knepley } else if (by == DM_BOUNDARY_MIRROR) { 5888865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(y - i - 1) + j; 58947c6ae99SBarry Smith } 590ac119b13SBarry Smith 59147c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 592ce00eea3SSatish Balay x_t = lx[n8 % m]; 59347c6ae99SBarry Smith /* y_t = ly[(n8/m)]; */ 59447c6ae99SBarry Smith s_t = bases[n8] + (i-1)*x_t; 5958865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 59647c6ae99SBarry Smith } 59747c6ae99SBarry Smith } 59847c6ae99SBarry Smith 5999566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(comm,dof,nn,idx,PETSC_USE_POINTER,&from)); 6009566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(global,from,local,to,>ol)); 6019566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)da,(PetscObject)gtol)); 6029566063dSJacob Faibussowitsch PetscCall(ISDestroy(&to)); 6039566063dSJacob Faibussowitsch PetscCall(ISDestroy(&from)); 60447c6ae99SBarry Smith 605d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_STAR) { 606ce00eea3SSatish Balay n0 = sn0; n2 = sn2; n6 = sn6; n8 = sn8; 607ce00eea3SSatish Balay } 608ce00eea3SSatish Balay 609288e7d53SBarry Smith if (((stencil_type == DMDA_STENCIL_STAR) || (bx && bx != DM_BOUNDARY_PERIODIC) || (by && by != DM_BOUNDARY_PERIODIC))) { 61047c6ae99SBarry Smith /* 61147c6ae99SBarry Smith Recompute the local to global mappings, this time keeping the 612ce00eea3SSatish Balay information about the cross corner processor numbers and any ghosted 613ce00eea3SSatish Balay but not periodic indices. 61447c6ae99SBarry Smith */ 61547c6ae99SBarry Smith nn = 0; 61647c6ae99SBarry Smith xbase = bases[rank]; 61747c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 61847c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 619ce00eea3SSatish Balay x_t = lx[n0 % m]; 62047c6ae99SBarry Smith y_t = ly[(n0/m)]; 62147c6ae99SBarry Smith s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x; 6228865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 623ce00eea3SSatish Balay } else if (xs-Xs > 0 && ys-Ys > 0) { 6248865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 62547c6ae99SBarry Smith } 62647c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 62747c6ae99SBarry Smith x_t = x; 62847c6ae99SBarry Smith y_t = ly[(n1/m)]; 62947c6ae99SBarry Smith s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t; 6308865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 631ce00eea3SSatish Balay } else if (ys-Ys > 0) { 632bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_MIRROR) { 6338865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(s_y - i + 1) + j; 634624904c4SBarry Smith } else { 6358865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 63647c6ae99SBarry Smith } 637624904c4SBarry Smith } 63847c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 639ce00eea3SSatish Balay x_t = lx[n2 % m]; 64047c6ae99SBarry Smith y_t = ly[(n2/m)]; 64147c6ae99SBarry Smith s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t; 6428865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 643ce00eea3SSatish Balay } else if (Xe-xe> 0 && ys-Ys > 0) { 6448865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 64547c6ae99SBarry Smith } 64647c6ae99SBarry Smith } 64747c6ae99SBarry Smith 64847c6ae99SBarry Smith for (i=0; i<y; i++) { 64947c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 650ce00eea3SSatish Balay x_t = lx[n3 % m]; 65147c6ae99SBarry Smith /* y_t = y; */ 65247c6ae99SBarry Smith s_t = bases[n3] + (i+1)*x_t - s_x; 6538865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 654ce00eea3SSatish Balay } else if (xs-Xs > 0) { 655bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_MIRROR) { 6568865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*i + s_x - j; 657624904c4SBarry Smith } else { 6588865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 65947c6ae99SBarry Smith } 660624904c4SBarry Smith } 66147c6ae99SBarry Smith 6628865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = xbase++; /* interior */ 66347c6ae99SBarry Smith 66447c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 665ce00eea3SSatish Balay x_t = lx[n5 % m]; 66647c6ae99SBarry Smith /* y_t = y; */ 66747c6ae99SBarry Smith s_t = bases[n5] + (i)*x_t; 6688865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 669ce00eea3SSatish Balay } else if (Xe-xe > 0) { 670bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_MIRROR) { 6718865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*(i + 1) - 2 - j; 672624904c4SBarry Smith } else { 6738865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 67447c6ae99SBarry Smith } 67547c6ae99SBarry Smith } 676624904c4SBarry Smith } 67747c6ae99SBarry Smith 67847c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 67947c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 680ce00eea3SSatish Balay x_t = lx[n6 % m]; 68147c6ae99SBarry Smith /* y_t = ly[(n6/m)]; */ 68247c6ae99SBarry Smith s_t = bases[n6] + (i)*x_t - s_x; 6838865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 684ce00eea3SSatish Balay } else if (xs-Xs > 0 && Ye-ye > 0) { 6858865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 68647c6ae99SBarry Smith } 68747c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 68847c6ae99SBarry Smith x_t = x; 68947c6ae99SBarry Smith /* y_t = ly[(n7/m)]; */ 69047c6ae99SBarry Smith s_t = bases[n7] + (i-1)*x_t; 6918865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 692ce00eea3SSatish Balay } else if (Ye-ye > 0) { 693bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_MIRROR) { 6948865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(y - i - 1) + j; 695624904c4SBarry Smith } else { 6968865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 69747c6ae99SBarry Smith } 698624904c4SBarry Smith } 69947c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 700ce00eea3SSatish Balay x_t = lx[n8 % m]; 70147c6ae99SBarry Smith /* y_t = ly[(n8/m)]; */ 70247c6ae99SBarry Smith s_t = bases[n8] + (i-1)*x_t; 7038865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 704ce00eea3SSatish Balay } else if (Xe-xe > 0 && Ye-ye > 0) { 7058865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 70647c6ae99SBarry Smith } 70747c6ae99SBarry Smith } 70847c6ae99SBarry Smith } 709ce00eea3SSatish Balay /* 710ce00eea3SSatish Balay Set the local to global ordering in the global vector, this allows use 711ce00eea3SSatish Balay of VecSetValuesLocal(). 712ce00eea3SSatish Balay */ 7139566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(comm,dof,nn,idx,PETSC_OWN_POINTER,&da->ltogmap)); 7149566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)da,(PetscObject)da->ltogmap)); 71547c6ae99SBarry Smith 7169566063dSJacob Faibussowitsch PetscCall(PetscFree2(bases,ldims)); 71747c6ae99SBarry Smith dd->m = m; dd->n = n; 718ce00eea3SSatish Balay /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */ 719ce00eea3SSatish Balay dd->xs = xs*dof; dd->xe = xe*dof; dd->ys = ys; dd->ye = ye; dd->zs = 0; dd->ze = 1; 720ce00eea3SSatish Balay dd->Xs = Xs*dof; dd->Xe = Xe*dof; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = 0; dd->Ze = 1; 72147c6ae99SBarry Smith 7229566063dSJacob Faibussowitsch PetscCall(VecDestroy(&local)); 7239566063dSJacob Faibussowitsch PetscCall(VecDestroy(&global)); 72447c6ae99SBarry Smith 72547c6ae99SBarry Smith dd->gtol = gtol; 72647c6ae99SBarry Smith dd->base = base; 7279a42bb27SBarry Smith da->ops->view = DMView_DA_2d; 7280298fd71SBarry Smith dd->ltol = NULL; 7290298fd71SBarry Smith dd->ao = NULL; 73047c6ae99SBarry Smith PetscFunctionReturn(0); 73147c6ae99SBarry Smith } 73247c6ae99SBarry Smith 73347c6ae99SBarry Smith /*@C 734aa219208SBarry Smith DMDACreate2d - Creates an object that will manage the communication of two-dimensional 73547c6ae99SBarry Smith regular array data that is distributed across some processors. 73647c6ae99SBarry Smith 737d083f849SBarry Smith Collective 73847c6ae99SBarry Smith 73947c6ae99SBarry Smith Input Parameters: 74047c6ae99SBarry Smith + comm - MPI communicator 7411321219cSEthan Coon . bx,by - type of ghost nodes the array have. 742bff4a2f0SMatthew G. Knepley Use one of DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_PERIODIC. 743aa219208SBarry Smith . stencil_type - stencil type. Use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR. 744897f7067SBarry Smith . M,N - global dimension in each direction of the array 74547c6ae99SBarry Smith . m,n - corresponding number of processors in each dimension 74647c6ae99SBarry Smith (or PETSC_DECIDE to have calculated) 74747c6ae99SBarry Smith . dof - number of degrees of freedom per node 74847c6ae99SBarry Smith . s - stencil width 74947c6ae99SBarry Smith - lx, ly - arrays containing the number of nodes in each cell along 7500298fd71SBarry Smith the x and y coordinates, or NULL. If non-null, these 75147c6ae99SBarry Smith must be of length as m and n, and the corresponding 75247c6ae99SBarry Smith m and n cannot be PETSC_DECIDE. The sum of the lx[] entries 75347c6ae99SBarry Smith must be M, and the sum of the ly[] entries must be N. 75447c6ae99SBarry Smith 75547c6ae99SBarry Smith Output Parameter: 75647c6ae99SBarry Smith . da - the resulting distributed array object 75747c6ae99SBarry Smith 75847c6ae99SBarry Smith Options Database Key: 759706a11cbSBarry Smith + -dm_view - Calls DMView() at the conclusion of DMDACreate2d() 760e43dc8daSBarry Smith . -da_grid_x <nx> - number of grid points in x direction 761e43dc8daSBarry Smith . -da_grid_y <ny> - number of grid points in y direction 76247c6ae99SBarry Smith . -da_processors_x <nx> - number of processors in x direction 76347c6ae99SBarry Smith . -da_processors_y <ny> - number of processors in y direction 764e0f5d30fSBarry Smith . -da_refine_x <rx> - refinement ratio in x direction 765e0f5d30fSBarry Smith . -da_refine_y <ry> - refinement ratio in y direction 766e43dc8daSBarry Smith - -da_refine <n> - refine the DMDA n times before creating 767e0f5d30fSBarry Smith 76847c6ae99SBarry Smith Level: beginner 76947c6ae99SBarry Smith 77047c6ae99SBarry Smith Notes: 771aa219208SBarry Smith The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the 772aa219208SBarry Smith standard 5-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes 77347c6ae99SBarry Smith the standard 9-pt stencil. 77447c6ae99SBarry Smith 775aa219208SBarry Smith The array data itself is NOT stored in the DMDA, it is stored in Vec objects; 776564755cdSBarry Smith The appropriate vector objects can be obtained with calls to DMCreateGlobalVector() 777564755cdSBarry Smith and DMCreateLocalVector() and calls to VecDuplicate() if more are needed. 77847c6ae99SBarry Smith 779897f7067SBarry Smith You must call DMSetUp() after this call before using this DM. 780897f7067SBarry Smith 781897f7067SBarry Smith If you wish to use the options database to change values in the DMDA call DMSetFromOptions() after this call 782897f7067SBarry Smith but before DMSetUp(). 783897f7067SBarry Smith 784db781477SPatrick Sanan .seealso: `DMDestroy()`, `DMView()`, `DMDACreate1d()`, `DMDACreate3d()`, `DMGlobalToLocalBegin()`, `DMDAGetRefinementFactor()`, 785db781477SPatrick Sanan `DMGlobalToLocalEnd()`, `DMLocalToGlobalBegin()`, `DMLocalToLocalBegin()`, `DMLocalToLocalEnd()`, `DMDASetRefinementFactor()`, 786db781477SPatrick Sanan `DMDAGetInfo()`, `DMCreateGlobalVector()`, `DMCreateLocalVector()`, `DMDACreateNaturalVector()`, `DMLoad()`, `DMDAGetOwnershipRanges()`, 787db781477SPatrick Sanan `DMStagCreate2d()` 78847c6ae99SBarry Smith 78947c6ae99SBarry Smith @*/ 790fe16a2e9SBarry Smith 791bff4a2f0SMatthew G. Knepley PetscErrorCode DMDACreate2d(MPI_Comm comm,DMBoundaryType bx,DMBoundaryType by,DMDAStencilType stencil_type, 7929a42bb27SBarry Smith PetscInt M,PetscInt N,PetscInt m,PetscInt n,PetscInt dof,PetscInt s,const PetscInt lx[],const PetscInt ly[],DM *da) 79347c6ae99SBarry Smith { 79447c6ae99SBarry Smith PetscFunctionBegin; 7959566063dSJacob Faibussowitsch PetscCall(DMDACreate(comm, da)); 7969566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*da, 2)); 7979566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(*da, M, N, 1)); 7989566063dSJacob Faibussowitsch PetscCall(DMDASetNumProcs(*da, m, n, PETSC_DECIDE)); 7999566063dSJacob Faibussowitsch PetscCall(DMDASetBoundaryType(*da, bx, by, DM_BOUNDARY_NONE)); 8009566063dSJacob Faibussowitsch PetscCall(DMDASetDof(*da, dof)); 8019566063dSJacob Faibussowitsch PetscCall(DMDASetStencilType(*da, stencil_type)); 8029566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(*da, s)); 8039566063dSJacob Faibussowitsch PetscCall(DMDASetOwnershipRanges(*da, lx, ly, NULL)); 80447c6ae99SBarry Smith PetscFunctionReturn(0); 80547c6ae99SBarry Smith } 806