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; 449566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Load Balance - Grid Points: Min %D avg %D max %D\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)); 519566063dSJacob Faibussowitsch PetscCall(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)); 529566063dSJacob Faibussowitsch PetscCall(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)); 539566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 549566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 558135c375SStefano Zampini } else if (format == PETSC_VIEWER_ASCII_GLVIS) { 569566063dSJacob Faibussowitsch PetscCall(DMView_DA_GLVis(da,viewer)); 573da9ae13SJed Brown } else { 589566063dSJacob Faibussowitsch PetscCall(DMView_DA_VTK(da,viewer)); 5947c6ae99SBarry Smith } 6047c6ae99SBarry Smith } else if (isdraw) { 6147c6ae99SBarry Smith PetscDraw draw; 6247c6ae99SBarry Smith double ymin = -1*dd->s-1,ymax = dd->N+dd->s; 6347c6ae99SBarry Smith double xmin = -1*dd->s-1,xmax = dd->M+dd->s; 6447c6ae99SBarry Smith double x,y; 658ea3bf28SBarry Smith PetscInt base; 668ea3bf28SBarry Smith const PetscInt *idx; 6747c6ae99SBarry Smith char node[10]; 6847c6ae99SBarry Smith PetscBool isnull; 6947c6ae99SBarry Smith 709566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(viewer,0,&draw)); 719566063dSJacob Faibussowitsch PetscCall(PetscDrawIsNull(draw,&isnull)); 725b399a63SLisandro Dalcin if (isnull) PetscFunctionReturn(0); 7347c6ae99SBarry Smith 749566063dSJacob Faibussowitsch PetscCall(PetscDrawCheckResizedWindow(draw)); 759566063dSJacob Faibussowitsch PetscCall(PetscDrawClear(draw)); 769566063dSJacob Faibussowitsch PetscCall(PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax)); 775b399a63SLisandro Dalcin 78*d0609cedSBarry Smith PetscDrawCollectiveBegin(draw); 7947c6ae99SBarry Smith /* first processor draw all node lines */ 80dd400576SPatrick Sanan if (rank == 0) { 8147c6ae99SBarry Smith ymin = 0.0; ymax = dd->N - 1; 8247c6ae99SBarry Smith for (xmin=0; xmin<dd->M; xmin++) { 839566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK)); 8447c6ae99SBarry Smith } 8547c6ae99SBarry Smith xmin = 0.0; xmax = dd->M - 1; 8647c6ae99SBarry Smith for (ymin=0; ymin<dd->N; ymin++) { 879566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK)); 8847c6ae99SBarry Smith } 8947c6ae99SBarry Smith } 90*d0609cedSBarry Smith PetscDrawCollectiveEnd(draw); 919566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 929566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw)); 9347c6ae99SBarry Smith 94*d0609cedSBarry Smith PetscDrawCollectiveBegin(draw); 9547c6ae99SBarry Smith /* draw my box */ 965b399a63SLisandro Dalcin xmin = dd->xs/dd->w; xmax =(dd->xe-1)/dd->w; ymin = dd->ys; ymax = dd->ye - 1; 979566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED)); 989566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED)); 999566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED)); 1009566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED)); 10147c6ae99SBarry Smith /* put in numbers */ 10247c6ae99SBarry Smith base = (dd->base)/dd->w; 10347c6ae99SBarry Smith for (y=ymin; y<=ymax; y++) { 10447c6ae99SBarry Smith for (x=xmin; x<=xmax; x++) { 1059566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(node,sizeof(node),"%d",(int)base++)); 1069566063dSJacob Faibussowitsch PetscCall(PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node)); 10747c6ae99SBarry Smith } 10847c6ae99SBarry Smith } 109*d0609cedSBarry Smith PetscDrawCollectiveEnd(draw); 1109566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 1119566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw)); 11247c6ae99SBarry Smith 113*d0609cedSBarry Smith PetscDrawCollectiveBegin(draw); 1145b399a63SLisandro Dalcin /* overlay ghost numbers, useful for error checking */ 1159566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(da->ltogmap,&idx)); 1165b399a63SLisandro Dalcin base = 0; xmin = dd->Xs; xmax = dd->Xe; ymin = dd->Ys; ymax = dd->Ye; 11747c6ae99SBarry Smith for (y=ymin; y<ymax; y++) { 11847c6ae99SBarry Smith for (x=xmin; x<xmax; x++) { 11947c6ae99SBarry Smith if ((base % dd->w) == 0) { 1209566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(node,sizeof(node),"%d",(int)(idx[base/dd->w]))); 1219566063dSJacob Faibussowitsch PetscCall(PetscDrawString(draw,x/dd->w,y,PETSC_DRAW_BLUE,node)); 12247c6ae99SBarry Smith } 12347c6ae99SBarry Smith base++; 12447c6ae99SBarry Smith } 12547c6ae99SBarry Smith } 1269566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(da->ltogmap,&idx)); 127*d0609cedSBarry Smith PetscDrawCollectiveEnd(draw); 1289566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 1299566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw)); 1309566063dSJacob Faibussowitsch PetscCall(PetscDrawSave(draw)); 131c9493c35SStefano Zampini } else if (isglvis) { 1329566063dSJacob Faibussowitsch PetscCall(DMView_DA_GLVis(da,viewer)); 1339a42bb27SBarry Smith } else if (isbinary) { 1349566063dSJacob Faibussowitsch PetscCall(DMView_DA_Binary(da,viewer)); 1359a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 1369a42bb27SBarry Smith } else if (ismatlab) { 1379566063dSJacob Faibussowitsch PetscCall(DMView_DA_Matlab(da,viewer)); 1389a42bb27SBarry Smith #endif 13911aeaf0aSBarry Smith } 14047c6ae99SBarry Smith PetscFunctionReturn(0); 14147c6ae99SBarry Smith } 14247c6ae99SBarry Smith 14347c6ae99SBarry Smith #if defined(new) 14447c6ae99SBarry Smith /* 145aa219208SBarry Smith DMDAGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix where local 146aa219208SBarry Smith function lives on a DMDA 14747c6ae99SBarry Smith 14847c6ae99SBarry Smith y ~= (F(u + ha) - F(u))/h, 14947c6ae99SBarry Smith where F = nonlinear function, as set by SNESSetFunction() 15047c6ae99SBarry Smith u = current iterate 15147c6ae99SBarry Smith h = difference interval 15247c6ae99SBarry Smith */ 153aa219208SBarry Smith PetscErrorCode DMDAGetDiagonal_MFFD(DM da,Vec U,Vec a) 15447c6ae99SBarry Smith { 15547c6ae99SBarry Smith PetscScalar h,*aa,*ww,v; 15647c6ae99SBarry Smith PetscReal epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON; 15747c6ae99SBarry Smith PetscInt gI,nI; 15847c6ae99SBarry Smith MatStencil stencil; 159aa219208SBarry Smith DMDALocalInfo info; 16047c6ae99SBarry Smith 16147c6ae99SBarry Smith PetscFunctionBegin; 1629566063dSJacob Faibussowitsch PetscCall((*ctx->func)(0,U,a,ctx->funcctx)); 1639566063dSJacob Faibussowitsch PetscCall((*ctx->funcisetbase)(U,ctx->funcctx)); 16447c6ae99SBarry Smith 1659566063dSJacob Faibussowitsch PetscCall(VecGetArray(U,&ww)); 1669566063dSJacob Faibussowitsch PetscCall(VecGetArray(a,&aa)); 16747c6ae99SBarry Smith 16847c6ae99SBarry Smith nI = 0; 16947c6ae99SBarry Smith h = ww[gI]; 17047c6ae99SBarry Smith if (h == 0.0) h = 1.0; 17147c6ae99SBarry Smith if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin; 17247c6ae99SBarry Smith else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin; 17347c6ae99SBarry Smith h *= epsilon; 17447c6ae99SBarry Smith 17547c6ae99SBarry Smith ww[gI] += h; 1769566063dSJacob Faibussowitsch PetscCall((*ctx->funci)(i,w,&v,ctx->funcctx)); 17747c6ae99SBarry Smith aa[nI] = (v - aa[nI])/h; 17847c6ae99SBarry Smith ww[gI] -= h; 17947c6ae99SBarry Smith nI++; 1808865f1eaSKarl Rupp 1819566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U,&ww)); 1829566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(a,&aa)); 18347c6ae99SBarry Smith PetscFunctionReturn(0); 18447c6ae99SBarry Smith } 18547c6ae99SBarry Smith #endif 18647c6ae99SBarry Smith 1877087cfbeSBarry Smith PetscErrorCode DMSetUp_DA_2D(DM da) 18847c6ae99SBarry Smith { 18947c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 19047c6ae99SBarry Smith const PetscInt M = dd->M; 19147c6ae99SBarry Smith const PetscInt N = dd->N; 19247c6ae99SBarry Smith PetscInt m = dd->m; 19347c6ae99SBarry Smith PetscInt n = dd->n; 19447c6ae99SBarry Smith const PetscInt dof = dd->w; 19547c6ae99SBarry Smith const PetscInt s = dd->s; 196bff4a2f0SMatthew G. Knepley DMBoundaryType bx = dd->bx; 197bff4a2f0SMatthew G. Knepley DMBoundaryType by = dd->by; 19819fd82e9SBarry Smith DMDAStencilType stencil_type = dd->stencil_type; 19947c6ae99SBarry Smith PetscInt *lx = dd->lx; 20047c6ae99SBarry Smith PetscInt *ly = dd->ly; 20147c6ae99SBarry Smith MPI_Comm comm; 20247c6ae99SBarry Smith PetscMPIInt rank,size; 203bd1fc5aeSBarry Smith PetscInt xs,xe,ys,ye,x,y,Xs,Xe,Ys,Ye,IXs,IXe,IYs,IYe; 2048ea3bf28SBarry Smith PetscInt up,down,left,right,i,n0,n1,n2,n3,n5,n6,n7,n8,*idx,nn; 205db87c5ecSEthan Coon PetscInt xbase,*bases,*ldims,j,x_t,y_t,s_t,base,count; 20647c6ae99SBarry Smith PetscInt s_x,s_y; /* s proportionalized to w */ 20747c6ae99SBarry Smith PetscInt sn0 = 0,sn2 = 0,sn6 = 0,sn8 = 0; 20847c6ae99SBarry Smith Vec local,global; 209bd1fc5aeSBarry Smith VecScatter gtol; 21045b6f7e9SBarry Smith IS to,from; 21147c6ae99SBarry Smith 21247c6ae99SBarry Smith PetscFunctionBegin; 2137a8be351SBarry 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"); 2149566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)da,&comm)); 2153855c12bSBarry Smith #if !defined(PETSC_USE_64BIT_INDICES) 2167a8be351SBarry Smith PetscCheck(((PetscInt64) M)*((PetscInt64) N)*((PetscInt64) dof) <= (PetscInt64) PETSC_MPI_INT_MAX,comm,PETSC_ERR_INT_OVERFLOW,"Mesh of %D by %D by %D (dof) is too large for 32 bit indices",M,N,dof); 2173855c12bSBarry Smith #endif 21847c6ae99SBarry Smith 2199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 2209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm,&rank)); 22147c6ae99SBarry Smith 2227d310018SBarry Smith dd->p = 1; 22347c6ae99SBarry Smith if (m != PETSC_DECIDE) { 2247a8be351SBarry Smith PetscCheck(m >= 1,comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m); 2257a8be351SBarry Smith else PetscCheck(m <= size,comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size); 22647c6ae99SBarry Smith } 22747c6ae99SBarry Smith if (n != PETSC_DECIDE) { 2287a8be351SBarry Smith PetscCheck(n >= 1,comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n); 2297a8be351SBarry Smith else PetscCheck(n <= size,comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size); 23047c6ae99SBarry Smith } 23147c6ae99SBarry Smith 23247c6ae99SBarry Smith if (m == PETSC_DECIDE || n == PETSC_DECIDE) { 23347c6ae99SBarry Smith if (n != PETSC_DECIDE) { 23447c6ae99SBarry Smith m = size/n; 23547c6ae99SBarry Smith } else if (m != PETSC_DECIDE) { 23647c6ae99SBarry Smith n = size/m; 23747c6ae99SBarry Smith } else { 23847c6ae99SBarry Smith /* try for squarish distribution */ 239369cc0aeSBarry Smith m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)N))); 24047c6ae99SBarry Smith if (!m) m = 1; 24147c6ae99SBarry Smith while (m > 0) { 24247c6ae99SBarry Smith n = size/m; 24347c6ae99SBarry Smith if (m*n == size) break; 24447c6ae99SBarry Smith m--; 24547c6ae99SBarry Smith } 24647c6ae99SBarry Smith if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;} 24747c6ae99SBarry Smith } 2487a8be351SBarry Smith PetscCheck(m*n == size,comm,PETSC_ERR_PLIB,"Unable to create partition, check the size of the communicator and input m and n "); 2497a8be351SBarry Smith } else PetscCheck(m*n == size,comm,PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition"); 25047c6ae99SBarry Smith 2517a8be351SBarry Smith PetscCheck(M >= m,comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m); 2527a8be351SBarry Smith PetscCheck(N >= n,comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n); 25347c6ae99SBarry Smith 25447c6ae99SBarry Smith /* 25547c6ae99SBarry Smith Determine locally owned region 25647c6ae99SBarry Smith xs is the first local node number, x is the number of local nodes 25747c6ae99SBarry Smith */ 25847c6ae99SBarry Smith if (!lx) { 2599566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m, &dd->lx)); 26047c6ae99SBarry Smith lx = dd->lx; 26147c6ae99SBarry Smith for (i=0; i<m; i++) { 26247c6ae99SBarry Smith lx[i] = M/m + ((M % m) > i); 26347c6ae99SBarry Smith } 26447c6ae99SBarry Smith } 26547c6ae99SBarry Smith x = lx[rank % m]; 26647c6ae99SBarry Smith xs = 0; 26747c6ae99SBarry Smith for (i=0; i<(rank % m); i++) { 26847c6ae99SBarry Smith xs += lx[i]; 26947c6ae99SBarry Smith } 27076bd3646SJed Brown if (PetscDefined(USE_DEBUG)) { 27147c6ae99SBarry Smith left = xs; 27247c6ae99SBarry Smith for (i=(rank % m); i<m; i++) { 27347c6ae99SBarry Smith left += lx[i]; 27447c6ae99SBarry Smith } 2757a8be351SBarry Smith PetscCheck(left == M,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Sum of lx across processors not equal to M: %D %D",left,M); 27676bd3646SJed Brown } 27747c6ae99SBarry Smith 27847c6ae99SBarry Smith /* 27947c6ae99SBarry Smith Determine locally owned region 28047c6ae99SBarry Smith ys is the first local node number, y is the number of local nodes 28147c6ae99SBarry Smith */ 28247c6ae99SBarry Smith if (!ly) { 2839566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &dd->ly)); 28447c6ae99SBarry Smith ly = dd->ly; 28547c6ae99SBarry Smith for (i=0; i<n; i++) { 28647c6ae99SBarry Smith ly[i] = N/n + ((N % n) > i); 28747c6ae99SBarry Smith } 28847c6ae99SBarry Smith } 28947c6ae99SBarry Smith y = ly[rank/m]; 29047c6ae99SBarry Smith ys = 0; 29147c6ae99SBarry Smith for (i=0; i<(rank/m); i++) { 29247c6ae99SBarry Smith ys += ly[i]; 29347c6ae99SBarry Smith } 29476bd3646SJed Brown if (PetscDefined(USE_DEBUG)) { 29547c6ae99SBarry Smith left = ys; 29647c6ae99SBarry Smith for (i=(rank/m); i<n; i++) { 29747c6ae99SBarry Smith left += ly[i]; 29847c6ae99SBarry Smith } 2997a8be351SBarry Smith PetscCheck(left == N,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Sum of ly across processors not equal to N: %D %D",left,N); 30076bd3646SJed Brown } 30147c6ae99SBarry Smith 302bcea557cSEthan Coon /* 303bcea557cSEthan Coon check if the scatter requires more than one process neighbor or wraps around 304bcea557cSEthan Coon the domain more than once 305bcea557cSEthan Coon */ 3067a8be351SBarry Smith PetscCheck((x >= s) || ((m <= 1) && (bx != DM_BOUNDARY_PERIODIC)),PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local x-width of domain x %D is smaller than stencil width s %D",x,s); 3077a8be351SBarry Smith PetscCheck((y >= s) || ((n <= 1) && (by != DM_BOUNDARY_PERIODIC)),PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local y-width of domain y %D is smaller than stencil width s %D",y,s); 30847c6ae99SBarry Smith xe = xs + x; 30947c6ae99SBarry Smith ye = ys + y; 31047c6ae99SBarry Smith 311ce00eea3SSatish Balay /* determine ghost region (Xs) and region scattered into (IXs) */ 312d9c9ebe5SPeter Brune if (xs-s > 0) { 313d9c9ebe5SPeter Brune Xs = xs - s; IXs = xs - s; 31488661749SPeter Brune } else { 31588661749SPeter Brune if (bx) { 31688661749SPeter Brune Xs = xs - s; 31788661749SPeter Brune } else { 31888661749SPeter Brune Xs = 0; 31988661749SPeter Brune } 32088661749SPeter Brune IXs = 0; 32188661749SPeter Brune } 322d9c9ebe5SPeter Brune if (xe+s <= M) { 323d9c9ebe5SPeter Brune Xe = xe + s; IXe = xe + s; 32488661749SPeter Brune } else { 32588661749SPeter Brune if (bx) { 326d9c9ebe5SPeter Brune Xs = xs - s; Xe = xe + s; 32788661749SPeter Brune } else { 32888661749SPeter Brune Xe = M; 32988661749SPeter Brune } 33088661749SPeter Brune IXe = M; 33188661749SPeter Brune } 33247c6ae99SBarry Smith 333bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC || bx == DM_BOUNDARY_MIRROR) { 334d9c9ebe5SPeter Brune IXs = xs - s; 335d9c9ebe5SPeter Brune IXe = xe + s; 336d9c9ebe5SPeter Brune Xs = xs - s; 337d9c9ebe5SPeter Brune Xe = xe + s; 33888661749SPeter Brune } 33947c6ae99SBarry Smith 340d9c9ebe5SPeter Brune if (ys-s > 0) { 341d9c9ebe5SPeter Brune Ys = ys - s; IYs = ys - s; 34288661749SPeter Brune } else { 34388661749SPeter Brune if (by) { 34488661749SPeter Brune Ys = ys - s; 34588661749SPeter Brune } else { 34688661749SPeter Brune Ys = 0; 34788661749SPeter Brune } 34888661749SPeter Brune IYs = 0; 34988661749SPeter Brune } 350d9c9ebe5SPeter Brune if (ye+s <= N) { 351d9c9ebe5SPeter Brune Ye = ye + s; IYe = ye + s; 35288661749SPeter Brune } else { 35388661749SPeter Brune if (by) { 35488661749SPeter Brune Ye = ye + s; 35588661749SPeter Brune } else { 35688661749SPeter Brune Ye = N; 35788661749SPeter Brune } 35888661749SPeter Brune IYe = N; 35988661749SPeter Brune } 36088661749SPeter Brune 361bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_PERIODIC || by == DM_BOUNDARY_MIRROR) { 362d9c9ebe5SPeter Brune IYs = ys - s; 363d9c9ebe5SPeter Brune IYe = ye + s; 364d9c9ebe5SPeter Brune Ys = ys - s; 365d9c9ebe5SPeter Brune Ye = ye + s; 36688661749SPeter Brune } 36788661749SPeter Brune 36888661749SPeter Brune /* stencil length in each direction */ 369d9c9ebe5SPeter Brune s_x = s; 370d9c9ebe5SPeter Brune s_y = s; 37147c6ae99SBarry Smith 37247c6ae99SBarry Smith /* determine starting point of each processor */ 37347c6ae99SBarry Smith nn = x*y; 3749566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(size+1,&bases,size,&ldims)); 3759566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm)); 37647c6ae99SBarry Smith bases[0] = 0; 37747c6ae99SBarry Smith for (i=1; i<=size; i++) { 37847c6ae99SBarry Smith bases[i] = ldims[i-1]; 37947c6ae99SBarry Smith } 38047c6ae99SBarry Smith for (i=1; i<=size; i++) { 38147c6ae99SBarry Smith bases[i] += bases[i-1]; 38247c6ae99SBarry Smith } 383ce00eea3SSatish Balay base = bases[rank]*dof; 38447c6ae99SBarry Smith 38547c6ae99SBarry Smith /* allocate the base parallel and sequential vectors */ 386ce00eea3SSatish Balay dd->Nlocal = x*y*dof; 3879566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,NULL,&global)); 388ce00eea3SSatish Balay dd->nlocal = (Xe-Xs)*(Ye-Ys)*dof; 3899566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,NULL,&local)); 39047c6ae99SBarry Smith 3914104a7a0SPatrick Sanan /* generate global to local vector scatter and local to global mapping*/ 39247c6ae99SBarry Smith 393ce00eea3SSatish Balay /* global to local must include ghost points within the domain, 394ce00eea3SSatish Balay but not ghost points outside the domain that aren't periodic */ 3959566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((IXe-IXs)*(IYe-IYs),&idx)); 396d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_BOX) { 397ce00eea3SSatish Balay left = IXs - Xs; right = left + (IXe-IXs); 398ce00eea3SSatish Balay down = IYs - Ys; up = down + (IYe-IYs); 399ce00eea3SSatish Balay count = 0; 400ce00eea3SSatish Balay for (i=down; i<up; i++) { 401ce00eea3SSatish Balay for (j=left; j<right; j++) { 402ce00eea3SSatish Balay idx[count++] = j + i*(Xe-Xs); 403ce00eea3SSatish Balay } 404ce00eea3SSatish Balay } 4059566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to)); 406ce00eea3SSatish Balay 40747c6ae99SBarry Smith } else { 40847c6ae99SBarry Smith /* must drop into cross shape region */ 40947c6ae99SBarry Smith /* ---------| 41047c6ae99SBarry Smith | top | 411ce00eea3SSatish Balay |--- ---| up 41247c6ae99SBarry Smith | middle | 41347c6ae99SBarry Smith | | 414ce00eea3SSatish Balay ---- ---- down 41547c6ae99SBarry Smith | bottom | 41647c6ae99SBarry Smith ----------- 41747c6ae99SBarry Smith Xs xs xe Xe */ 418ce00eea3SSatish Balay left = xs - Xs; right = left + x; 419ce00eea3SSatish Balay down = ys - Ys; up = down + y; 42047c6ae99SBarry Smith count = 0; 421ce00eea3SSatish Balay /* bottom */ 422ce00eea3SSatish Balay for (i=(IYs-Ys); i<down; i++) { 423ce00eea3SSatish Balay for (j=left; j<right; j++) { 424ce00eea3SSatish Balay idx[count++] = j + i*(Xe-Xs); 42547c6ae99SBarry Smith } 42647c6ae99SBarry Smith } 42747c6ae99SBarry Smith /* middle */ 42847c6ae99SBarry Smith for (i=down; i<up; i++) { 429ce00eea3SSatish Balay for (j=(IXs-Xs); j<(IXe-Xs); j++) { 430ce00eea3SSatish Balay idx[count++] = j + i*(Xe-Xs); 43147c6ae99SBarry Smith } 43247c6ae99SBarry Smith } 43347c6ae99SBarry Smith /* top */ 434ce00eea3SSatish Balay for (i=up; i<up+IYe-ye; i++) { 435ce00eea3SSatish Balay for (j=left; j<right; j++) { 436ce00eea3SSatish Balay idx[count++] = j + i*(Xe-Xs); 43747c6ae99SBarry Smith } 43847c6ae99SBarry Smith } 4399566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to)); 44047c6ae99SBarry Smith } 44147c6ae99SBarry Smith 44247c6ae99SBarry Smith /* determine who lies on each side of us stored in n6 n7 n8 44347c6ae99SBarry Smith n3 n5 44447c6ae99SBarry Smith n0 n1 n2 44547c6ae99SBarry Smith */ 44647c6ae99SBarry Smith 44747c6ae99SBarry Smith /* Assume the Non-Periodic Case */ 44847c6ae99SBarry Smith n1 = rank - m; 44947c6ae99SBarry Smith if (rank % m) { 45047c6ae99SBarry Smith n0 = n1 - 1; 45147c6ae99SBarry Smith } else { 45247c6ae99SBarry Smith n0 = -1; 45347c6ae99SBarry Smith } 45447c6ae99SBarry Smith if ((rank+1) % m) { 45547c6ae99SBarry Smith n2 = n1 + 1; 45647c6ae99SBarry Smith n5 = rank + 1; 45747c6ae99SBarry Smith n8 = rank + m + 1; if (n8 >= m*n) n8 = -1; 45847c6ae99SBarry Smith } else { 45947c6ae99SBarry Smith n2 = -1; n5 = -1; n8 = -1; 46047c6ae99SBarry Smith } 46147c6ae99SBarry Smith if (rank % m) { 46247c6ae99SBarry Smith n3 = rank - 1; 46347c6ae99SBarry Smith n6 = n3 + m; if (n6 >= m*n) n6 = -1; 46447c6ae99SBarry Smith } else { 46547c6ae99SBarry Smith n3 = -1; n6 = -1; 46647c6ae99SBarry Smith } 46747c6ae99SBarry Smith n7 = rank + m; if (n7 >= m*n) n7 = -1; 46847c6ae99SBarry Smith 469bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC && by == DM_BOUNDARY_PERIODIC) { 47047c6ae99SBarry Smith /* Modify for Periodic Cases */ 47147c6ae99SBarry Smith /* Handle all four corners */ 47247c6ae99SBarry Smith if ((n6 < 0) && (n7 < 0) && (n3 < 0)) n6 = m-1; 47347c6ae99SBarry Smith if ((n8 < 0) && (n7 < 0) && (n5 < 0)) n8 = 0; 47447c6ae99SBarry Smith if ((n2 < 0) && (n5 < 0) && (n1 < 0)) n2 = size-m; 47547c6ae99SBarry Smith if ((n0 < 0) && (n3 < 0) && (n1 < 0)) n0 = size-1; 47647c6ae99SBarry Smith 47747c6ae99SBarry Smith /* Handle Top and Bottom Sides */ 47847c6ae99SBarry Smith if (n1 < 0) n1 = rank + m * (n-1); 47947c6ae99SBarry Smith if (n7 < 0) n7 = rank - m * (n-1); 48047c6ae99SBarry Smith if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1; 48147c6ae99SBarry Smith if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1; 48247c6ae99SBarry Smith if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1; 48347c6ae99SBarry Smith if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1; 48447c6ae99SBarry Smith 48547c6ae99SBarry Smith /* Handle Left and Right Sides */ 48647c6ae99SBarry Smith if (n3 < 0) n3 = rank + (m-1); 48747c6ae99SBarry Smith if (n5 < 0) n5 = rank - (m-1); 48847c6ae99SBarry Smith if ((n1 >= 0) && (n0 < 0)) n0 = rank-1; 48947c6ae99SBarry Smith if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1; 49047c6ae99SBarry Smith if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1; 49147c6ae99SBarry Smith if ((n7 >= 0) && (n8 < 0)) n8 = rank+1; 492bff4a2f0SMatthew G. Knepley } else if (by == DM_BOUNDARY_PERIODIC) { /* Handle Top and Bottom Sides */ 493ce00eea3SSatish Balay if (n1 < 0) n1 = rank + m * (n-1); 494ce00eea3SSatish Balay if (n7 < 0) n7 = rank - m * (n-1); 495ce00eea3SSatish Balay if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1; 496ce00eea3SSatish Balay if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1; 497ce00eea3SSatish Balay if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1; 498ce00eea3SSatish Balay if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1; 499bff4a2f0SMatthew G. Knepley } else if (bx == DM_BOUNDARY_PERIODIC) { /* Handle Left and Right Sides */ 500ce00eea3SSatish Balay if (n3 < 0) n3 = rank + (m-1); 501ce00eea3SSatish Balay if (n5 < 0) n5 = rank - (m-1); 502ce00eea3SSatish Balay if ((n1 >= 0) && (n0 < 0)) n0 = rank-1; 503ce00eea3SSatish Balay if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1; 504ce00eea3SSatish Balay if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1; 505ce00eea3SSatish Balay if ((n7 >= 0) && (n8 < 0)) n8 = rank+1; 50647c6ae99SBarry Smith } 507ce00eea3SSatish Balay 5089566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(9,&dd->neighbors)); 5098865f1eaSKarl Rupp 51047c6ae99SBarry Smith dd->neighbors[0] = n0; 51147c6ae99SBarry Smith dd->neighbors[1] = n1; 51247c6ae99SBarry Smith dd->neighbors[2] = n2; 51347c6ae99SBarry Smith dd->neighbors[3] = n3; 51447c6ae99SBarry Smith dd->neighbors[4] = rank; 51547c6ae99SBarry Smith dd->neighbors[5] = n5; 51647c6ae99SBarry Smith dd->neighbors[6] = n6; 51747c6ae99SBarry Smith dd->neighbors[7] = n7; 51847c6ae99SBarry Smith dd->neighbors[8] = n8; 51947c6ae99SBarry Smith 520d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_STAR) { 52147c6ae99SBarry Smith /* save corner processor numbers */ 52247c6ae99SBarry Smith sn0 = n0; sn2 = n2; sn6 = n6; sn8 = n8; 52347c6ae99SBarry Smith n0 = n2 = n6 = n8 = -1; 52447c6ae99SBarry Smith } 52547c6ae99SBarry Smith 5269566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((Xe-Xs)*(Ye-Ys),&idx)); 52747c6ae99SBarry Smith 528ce00eea3SSatish Balay nn = 0; 52947c6ae99SBarry Smith xbase = bases[rank]; 53047c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 53147c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 532ce00eea3SSatish Balay x_t = lx[n0 % m]; 53347c6ae99SBarry Smith y_t = ly[(n0/m)]; 53447c6ae99SBarry Smith s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x; 5358865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 53647c6ae99SBarry Smith } 537ac119b13SBarry Smith 53847c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 53947c6ae99SBarry Smith x_t = x; 54047c6ae99SBarry Smith y_t = ly[(n1/m)]; 54147c6ae99SBarry Smith s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t; 5428865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 543bff4a2f0SMatthew G. Knepley } else if (by == DM_BOUNDARY_MIRROR) { 5448865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(s_y - i + 1) + j; 54547c6ae99SBarry Smith } 546ac119b13SBarry Smith 54747c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 548ce00eea3SSatish Balay x_t = lx[n2 % m]; 54947c6ae99SBarry Smith y_t = ly[(n2/m)]; 55047c6ae99SBarry Smith s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t; 5518865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 55247c6ae99SBarry Smith } 55347c6ae99SBarry Smith } 55447c6ae99SBarry Smith 55547c6ae99SBarry Smith for (i=0; i<y; i++) { 55647c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 557ce00eea3SSatish Balay x_t = lx[n3 % m]; 55847c6ae99SBarry Smith /* y_t = y; */ 55947c6ae99SBarry Smith s_t = bases[n3] + (i+1)*x_t - s_x; 5608865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 561bff4a2f0SMatthew G. Knepley } else if (bx == DM_BOUNDARY_MIRROR) { 5628865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*i + s_x - j; 56347c6ae99SBarry Smith } 56447c6ae99SBarry Smith 5658865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = xbase++; /* interior */ 56647c6ae99SBarry Smith 56747c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 568ce00eea3SSatish Balay x_t = lx[n5 % m]; 56947c6ae99SBarry Smith /* y_t = y; */ 57047c6ae99SBarry Smith s_t = bases[n5] + (i)*x_t; 5718865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 572bff4a2f0SMatthew G. Knepley } else if (bx == DM_BOUNDARY_MIRROR) { 5738865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*(i + 1) - 2 - j; 57447c6ae99SBarry Smith } 57547c6ae99SBarry Smith } 57647c6ae99SBarry Smith 57747c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 57847c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 579ce00eea3SSatish Balay x_t = lx[n6 % m]; 58047c6ae99SBarry Smith /* y_t = ly[(n6/m)]; */ 58147c6ae99SBarry Smith s_t = bases[n6] + (i)*x_t - s_x; 5828865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 58347c6ae99SBarry Smith } 584ac119b13SBarry Smith 58547c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 58647c6ae99SBarry Smith x_t = x; 58747c6ae99SBarry Smith /* y_t = ly[(n7/m)]; */ 58847c6ae99SBarry Smith s_t = bases[n7] + (i-1)*x_t; 5898865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 590bff4a2f0SMatthew G. Knepley } else if (by == DM_BOUNDARY_MIRROR) { 5918865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(y - i - 1) + j; 59247c6ae99SBarry Smith } 593ac119b13SBarry Smith 59447c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 595ce00eea3SSatish Balay x_t = lx[n8 % m]; 59647c6ae99SBarry Smith /* y_t = ly[(n8/m)]; */ 59747c6ae99SBarry Smith s_t = bases[n8] + (i-1)*x_t; 5988865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 59947c6ae99SBarry Smith } 60047c6ae99SBarry Smith } 60147c6ae99SBarry Smith 6029566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(comm,dof,nn,idx,PETSC_USE_POINTER,&from)); 6039566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(global,from,local,to,>ol)); 6049566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)da,(PetscObject)gtol)); 6059566063dSJacob Faibussowitsch PetscCall(ISDestroy(&to)); 6069566063dSJacob Faibussowitsch PetscCall(ISDestroy(&from)); 60747c6ae99SBarry Smith 608d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_STAR) { 609ce00eea3SSatish Balay n0 = sn0; n2 = sn2; n6 = sn6; n8 = sn8; 610ce00eea3SSatish Balay } 611ce00eea3SSatish Balay 612288e7d53SBarry Smith if (((stencil_type == DMDA_STENCIL_STAR) || (bx && bx != DM_BOUNDARY_PERIODIC) || (by && by != DM_BOUNDARY_PERIODIC))) { 61347c6ae99SBarry Smith /* 61447c6ae99SBarry Smith Recompute the local to global mappings, this time keeping the 615ce00eea3SSatish Balay information about the cross corner processor numbers and any ghosted 616ce00eea3SSatish Balay but not periodic indices. 61747c6ae99SBarry Smith */ 61847c6ae99SBarry Smith nn = 0; 61947c6ae99SBarry Smith xbase = bases[rank]; 62047c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 62147c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 622ce00eea3SSatish Balay x_t = lx[n0 % m]; 62347c6ae99SBarry Smith y_t = ly[(n0/m)]; 62447c6ae99SBarry Smith s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x; 6258865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 626ce00eea3SSatish Balay } else if (xs-Xs > 0 && ys-Ys > 0) { 6278865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 62847c6ae99SBarry Smith } 62947c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 63047c6ae99SBarry Smith x_t = x; 63147c6ae99SBarry Smith y_t = ly[(n1/m)]; 63247c6ae99SBarry Smith s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t; 6338865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 634ce00eea3SSatish Balay } else if (ys-Ys > 0) { 635bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_MIRROR) { 6368865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(s_y - i + 1) + j; 637624904c4SBarry Smith } else { 6388865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 63947c6ae99SBarry Smith } 640624904c4SBarry Smith } 64147c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 642ce00eea3SSatish Balay x_t = lx[n2 % m]; 64347c6ae99SBarry Smith y_t = ly[(n2/m)]; 64447c6ae99SBarry Smith s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t; 6458865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 646ce00eea3SSatish Balay } else if (Xe-xe> 0 && ys-Ys > 0) { 6478865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 64847c6ae99SBarry Smith } 64947c6ae99SBarry Smith } 65047c6ae99SBarry Smith 65147c6ae99SBarry Smith for (i=0; i<y; i++) { 65247c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 653ce00eea3SSatish Balay x_t = lx[n3 % m]; 65447c6ae99SBarry Smith /* y_t = y; */ 65547c6ae99SBarry Smith s_t = bases[n3] + (i+1)*x_t - s_x; 6568865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 657ce00eea3SSatish Balay } else if (xs-Xs > 0) { 658bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_MIRROR) { 6598865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*i + s_x - j; 660624904c4SBarry Smith } else { 6618865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 66247c6ae99SBarry Smith } 663624904c4SBarry Smith } 66447c6ae99SBarry Smith 6658865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = xbase++; /* interior */ 66647c6ae99SBarry Smith 66747c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 668ce00eea3SSatish Balay x_t = lx[n5 % m]; 66947c6ae99SBarry Smith /* y_t = y; */ 67047c6ae99SBarry Smith s_t = bases[n5] + (i)*x_t; 6718865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 672ce00eea3SSatish Balay } else if (Xe-xe > 0) { 673bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_MIRROR) { 6748865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*(i + 1) - 2 - j; 675624904c4SBarry Smith } else { 6768865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 67747c6ae99SBarry Smith } 67847c6ae99SBarry Smith } 679624904c4SBarry Smith } 68047c6ae99SBarry Smith 68147c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 68247c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 683ce00eea3SSatish Balay x_t = lx[n6 % m]; 68447c6ae99SBarry Smith /* y_t = ly[(n6/m)]; */ 68547c6ae99SBarry Smith s_t = bases[n6] + (i)*x_t - s_x; 6868865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 687ce00eea3SSatish Balay } else if (xs-Xs > 0 && Ye-ye > 0) { 6888865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 68947c6ae99SBarry Smith } 69047c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 69147c6ae99SBarry Smith x_t = x; 69247c6ae99SBarry Smith /* y_t = ly[(n7/m)]; */ 69347c6ae99SBarry Smith s_t = bases[n7] + (i-1)*x_t; 6948865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 695ce00eea3SSatish Balay } else if (Ye-ye > 0) { 696bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_MIRROR) { 6978865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(y - i - 1) + j; 698624904c4SBarry Smith } else { 6998865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 70047c6ae99SBarry Smith } 701624904c4SBarry Smith } 70247c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 703ce00eea3SSatish Balay x_t = lx[n8 % m]; 70447c6ae99SBarry Smith /* y_t = ly[(n8/m)]; */ 70547c6ae99SBarry Smith s_t = bases[n8] + (i-1)*x_t; 7068865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 707ce00eea3SSatish Balay } else if (Xe-xe > 0 && Ye-ye > 0) { 7088865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 70947c6ae99SBarry Smith } 71047c6ae99SBarry Smith } 71147c6ae99SBarry Smith } 712ce00eea3SSatish Balay /* 713ce00eea3SSatish Balay Set the local to global ordering in the global vector, this allows use 714ce00eea3SSatish Balay of VecSetValuesLocal(). 715ce00eea3SSatish Balay */ 7169566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(comm,dof,nn,idx,PETSC_OWN_POINTER,&da->ltogmap)); 7179566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)da,(PetscObject)da->ltogmap)); 71847c6ae99SBarry Smith 7199566063dSJacob Faibussowitsch PetscCall(PetscFree2(bases,ldims)); 72047c6ae99SBarry Smith dd->m = m; dd->n = n; 721ce00eea3SSatish Balay /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */ 722ce00eea3SSatish Balay dd->xs = xs*dof; dd->xe = xe*dof; dd->ys = ys; dd->ye = ye; dd->zs = 0; dd->ze = 1; 723ce00eea3SSatish Balay dd->Xs = Xs*dof; dd->Xe = Xe*dof; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = 0; dd->Ze = 1; 72447c6ae99SBarry Smith 7259566063dSJacob Faibussowitsch PetscCall(VecDestroy(&local)); 7269566063dSJacob Faibussowitsch PetscCall(VecDestroy(&global)); 72747c6ae99SBarry Smith 72847c6ae99SBarry Smith dd->gtol = gtol; 72947c6ae99SBarry Smith dd->base = base; 7309a42bb27SBarry Smith da->ops->view = DMView_DA_2d; 7310298fd71SBarry Smith dd->ltol = NULL; 7320298fd71SBarry Smith dd->ao = NULL; 73347c6ae99SBarry Smith PetscFunctionReturn(0); 73447c6ae99SBarry Smith } 73547c6ae99SBarry Smith 73647c6ae99SBarry Smith /*@C 737aa219208SBarry Smith DMDACreate2d - Creates an object that will manage the communication of two-dimensional 73847c6ae99SBarry Smith regular array data that is distributed across some processors. 73947c6ae99SBarry Smith 740d083f849SBarry Smith Collective 74147c6ae99SBarry Smith 74247c6ae99SBarry Smith Input Parameters: 74347c6ae99SBarry Smith + comm - MPI communicator 7441321219cSEthan Coon . bx,by - type of ghost nodes the array have. 745bff4a2f0SMatthew G. Knepley Use one of DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_PERIODIC. 746aa219208SBarry Smith . stencil_type - stencil type. Use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR. 747897f7067SBarry Smith . M,N - global dimension in each direction of the array 74847c6ae99SBarry Smith . m,n - corresponding number of processors in each dimension 74947c6ae99SBarry Smith (or PETSC_DECIDE to have calculated) 75047c6ae99SBarry Smith . dof - number of degrees of freedom per node 75147c6ae99SBarry Smith . s - stencil width 75247c6ae99SBarry Smith - lx, ly - arrays containing the number of nodes in each cell along 7530298fd71SBarry Smith the x and y coordinates, or NULL. If non-null, these 75447c6ae99SBarry Smith must be of length as m and n, and the corresponding 75547c6ae99SBarry Smith m and n cannot be PETSC_DECIDE. The sum of the lx[] entries 75647c6ae99SBarry Smith must be M, and the sum of the ly[] entries must be N. 75747c6ae99SBarry Smith 75847c6ae99SBarry Smith Output Parameter: 75947c6ae99SBarry Smith . da - the resulting distributed array object 76047c6ae99SBarry Smith 76147c6ae99SBarry Smith Options Database Key: 762706a11cbSBarry Smith + -dm_view - Calls DMView() at the conclusion of DMDACreate2d() 763e43dc8daSBarry Smith . -da_grid_x <nx> - number of grid points in x direction 764e43dc8daSBarry Smith . -da_grid_y <ny> - number of grid points in y direction 76547c6ae99SBarry Smith . -da_processors_x <nx> - number of processors in x direction 76647c6ae99SBarry Smith . -da_processors_y <ny> - number of processors in y direction 767e0f5d30fSBarry Smith . -da_refine_x <rx> - refinement ratio in x direction 768e0f5d30fSBarry Smith . -da_refine_y <ry> - refinement ratio in y direction 769e43dc8daSBarry Smith - -da_refine <n> - refine the DMDA n times before creating 770e0f5d30fSBarry Smith 77147c6ae99SBarry Smith Level: beginner 77247c6ae99SBarry Smith 77347c6ae99SBarry Smith Notes: 774aa219208SBarry Smith The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the 775aa219208SBarry Smith standard 5-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes 77647c6ae99SBarry Smith the standard 9-pt stencil. 77747c6ae99SBarry Smith 778aa219208SBarry Smith The array data itself is NOT stored in the DMDA, it is stored in Vec objects; 779564755cdSBarry Smith The appropriate vector objects can be obtained with calls to DMCreateGlobalVector() 780564755cdSBarry Smith and DMCreateLocalVector() and calls to VecDuplicate() if more are needed. 78147c6ae99SBarry Smith 782897f7067SBarry Smith You must call DMSetUp() after this call before using this DM. 783897f7067SBarry Smith 784897f7067SBarry Smith If you wish to use the options database to change values in the DMDA call DMSetFromOptions() after this call 785897f7067SBarry Smith but before DMSetUp(). 786897f7067SBarry Smith 787aa219208SBarry Smith .seealso: DMDestroy(), DMView(), DMDACreate1d(), DMDACreate3d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(), 78899f0b487SRichard Tran Mills DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMLocalToLocalBegin(), DMLocalToLocalEnd(), DMDASetRefinementFactor(), 7893bd220d7SPatrick Sanan DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges(), 7903bd220d7SPatrick Sanan DMStagCreate2d() 79147c6ae99SBarry Smith 79247c6ae99SBarry Smith @*/ 793fe16a2e9SBarry Smith 794bff4a2f0SMatthew G. Knepley PetscErrorCode DMDACreate2d(MPI_Comm comm,DMBoundaryType bx,DMBoundaryType by,DMDAStencilType stencil_type, 7959a42bb27SBarry Smith PetscInt M,PetscInt N,PetscInt m,PetscInt n,PetscInt dof,PetscInt s,const PetscInt lx[],const PetscInt ly[],DM *da) 79647c6ae99SBarry Smith { 79747c6ae99SBarry Smith PetscFunctionBegin; 7989566063dSJacob Faibussowitsch PetscCall(DMDACreate(comm, da)); 7999566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*da, 2)); 8009566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(*da, M, N, 1)); 8019566063dSJacob Faibussowitsch PetscCall(DMDASetNumProcs(*da, m, n, PETSC_DECIDE)); 8029566063dSJacob Faibussowitsch PetscCall(DMDASetBoundaryType(*da, bx, by, DM_BOUNDARY_NONE)); 8039566063dSJacob Faibussowitsch PetscCall(DMDASetDof(*da, dof)); 8049566063dSJacob Faibussowitsch PetscCall(DMDASetStencilType(*da, stencil_type)); 8059566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(*da, s)); 8069566063dSJacob Faibussowitsch PetscCall(DMDASetOwnershipRanges(*da, lx, ly, NULL)); 80747c6ae99SBarry Smith PetscFunctionReturn(0); 80847c6ae99SBarry Smith } 809