1 2 /* 3 Code for manipulating distributed regular 1d arrays in parallel. 4 This file was created by Peter Mell 6/30/95 5 */ 6 7 #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 8 9 #include <petscdraw.h> 10 #undef __FUNCT__ 11 #define __FUNCT__ "DMView_DA_1d" 12 PetscErrorCode DMView_DA_1d(DM da,PetscViewer viewer) 13 { 14 PetscErrorCode ierr; 15 PetscMPIInt rank; 16 PetscBool iascii,isdraw,isbinary; 17 DM_DA *dd = (DM_DA*)da->data; 18 #if defined(PETSC_HAVE_MATLAB_ENGINE) 19 PetscBool ismatlab; 20 #endif 21 22 PetscFunctionBegin; 23 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)da),&rank);CHKERRQ(ierr); 24 25 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 26 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 27 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 28 #if defined(PETSC_HAVE_MATLAB_ENGINE) 29 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERMATLAB,&ismatlab);CHKERRQ(ierr); 30 #endif 31 if (iascii) { 32 PetscViewerFormat format; 33 34 ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr); 35 if (format != PETSC_VIEWER_ASCII_VTK && format != PETSC_VIEWER_ASCII_VTK_CELL) { 36 DMDALocalInfo info; 37 ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr); 38 ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 39 ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Processor [%d] M %D m %D w %D s %D\n",rank,dd->M,dd->m,dd->w,dd->s);CHKERRQ(ierr); 40 ierr = PetscViewerASCIISynchronizedPrintf(viewer,"X range of indices: %D %D\n",info.xs,info.xs+info.xm);CHKERRQ(ierr); 41 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 42 ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 43 } else { 44 ierr = DMView_DA_VTK(da, viewer);CHKERRQ(ierr); 45 } 46 } else if (isdraw) { 47 PetscDraw draw; 48 double ymin = -1,ymax = 1,xmin = -1,xmax = dd->M,x; 49 PetscInt base; 50 char node[10]; 51 PetscBool isnull; 52 53 ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 54 ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 55 if (isnull) PetscFunctionReturn(0); 56 57 ierr = PetscDrawCheckResizedWindow(draw);CHKERRQ(ierr); 58 ierr = PetscDrawClear(draw);CHKERRQ(ierr); 59 60 ierr = PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);CHKERRQ(ierr); 61 62 ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr); 63 /* first processor draws all node lines */ 64 if (!rank) { 65 PetscInt xmin_tmp; 66 ymin = 0.0; ymax = 0.3; 67 for (xmin_tmp=0; xmin_tmp < dd->M; xmin_tmp++) { 68 ierr = PetscDrawLine(draw,(double)xmin_tmp,ymin,(double)xmin_tmp,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr); 69 } 70 xmin = 0.0; xmax = dd->M - 1; 71 ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);CHKERRQ(ierr); 72 ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr); 73 } 74 ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr); 75 76 ierr = PetscDrawFlush(draw);CHKERRQ(ierr); 77 ierr = PetscDrawPause(draw);CHKERRQ(ierr); 78 79 ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr); 80 /* draw my box */ 81 ymin = 0; ymax = 0.3; xmin = dd->xs / dd->w; xmax = (dd->xe / dd->w) - 1; 82 ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);CHKERRQ(ierr); 83 ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 84 ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 85 ierr = PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 86 /* Put in index numbers */ 87 base = dd->base / dd->w; 88 for (x=xmin; x<=xmax; x++) { 89 sprintf(node,"%d",(int)base++); 90 ierr = PetscDrawString(draw,x,ymin,PETSC_DRAW_RED,node);CHKERRQ(ierr); 91 } 92 ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr); 93 94 ierr = PetscDrawFlush(draw);CHKERRQ(ierr); 95 ierr = PetscDrawPause(draw);CHKERRQ(ierr); 96 } else if (isbinary) { 97 ierr = DMView_DA_Binary(da,viewer);CHKERRQ(ierr); 98 #if defined(PETSC_HAVE_MATLAB_ENGINE) 99 } else if (ismatlab) { 100 ierr = DMView_DA_Matlab(da,viewer);CHKERRQ(ierr); 101 #endif 102 } 103 PetscFunctionReturn(0); 104 } 105 106 107 #undef __FUNCT__ 108 #define __FUNCT__ "DMSetUp_DA_1D" 109 PetscErrorCode DMSetUp_DA_1D(DM da) 110 { 111 DM_DA *dd = (DM_DA*)da->data; 112 const PetscInt M = dd->M; 113 const PetscInt dof = dd->w; 114 const PetscInt s = dd->s; 115 const PetscInt sDist = s; /* stencil distance in points */ 116 const PetscInt *lx = dd->lx; 117 DMBoundaryType bx = dd->bx; 118 MPI_Comm comm; 119 Vec local, global; 120 VecScatter gtol; 121 IS to, from; 122 PetscBool flg1 = PETSC_FALSE, flg2 = PETSC_FALSE; 123 PetscMPIInt rank, size; 124 PetscInt i,*idx,nn,left,xs,xe,x,Xs,Xe,start,m,IXs,IXe; 125 PetscErrorCode ierr; 126 127 PetscFunctionBegin; 128 ierr = PetscObjectGetComm((PetscObject) da, &comm);CHKERRQ(ierr); 129 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 130 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 131 132 dd->p = 1; 133 dd->n = 1; 134 dd->m = size; 135 m = dd->m; 136 137 if (s > 0) { 138 /* if not communicating data then should be ok to have nothing on some processes */ 139 if (M < m) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"More processes than data points! %D %D",m,M); 140 if ((M-1) < s && size > 1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array is too small for stencil! %D %D",M-1,s); 141 } 142 143 /* 144 Determine locally owned region 145 xs is the first local node number, x is the number of local nodes 146 */ 147 if (!lx) { 148 ierr = PetscMalloc1(m, &dd->lx);CHKERRQ(ierr); 149 ierr = PetscOptionsGetBool(((PetscObject)da)->options,((PetscObject)da)->prefix,"-da_partition_blockcomm",&flg1,NULL);CHKERRQ(ierr); 150 ierr = PetscOptionsGetBool(((PetscObject)da)->options,((PetscObject)da)->prefix,"-da_partition_nodes_at_end",&flg2,NULL);CHKERRQ(ierr); 151 if (flg1) { /* Block Comm type Distribution */ 152 xs = rank*M/m; 153 x = (rank + 1)*M/m - xs; 154 } else if (flg2) { /* The odd nodes are evenly distributed across last nodes */ 155 x = (M + rank)/m; 156 if (M/m == x) xs = rank*x; 157 else xs = rank*(x-1) + (M+rank)%(x*m); 158 } else { /* The odd nodes are evenly distributed across the first k nodes */ 159 /* Regular PETSc Distribution */ 160 x = M/m + ((M % m) > rank); 161 if (rank >= (M % m)) xs = (rank * (PetscInt)(M/m) + M % m); 162 else xs = rank * (PetscInt)(M/m) + rank; 163 } 164 ierr = MPI_Allgather(&xs,1,MPIU_INT,dd->lx,1,MPIU_INT,comm);CHKERRQ(ierr); 165 for (i=0; i<m-1; i++) dd->lx[i] = dd->lx[i+1] - dd->lx[i]; 166 dd->lx[m-1] = M - dd->lx[m-1]; 167 } else { 168 x = lx[rank]; 169 xs = 0; 170 for (i=0; i<rank; i++) xs += lx[i]; 171 /* verify that data user provided is consistent */ 172 left = xs; 173 for (i=rank; i<size; i++) left += lx[i]; 174 if (left != M) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Sum of lx across processors not equal to M %D %D",left,M); 175 } 176 177 /* 178 check if the scatter requires more than one process neighbor or wraps around 179 the domain more than once 180 */ 181 if ((x < s) & ((M > 1) | (bx == DM_BOUNDARY_PERIODIC))) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local x-width of domain x %D is smaller than stencil width s %D",x,s); 182 183 xe = xs + x; 184 185 /* determine ghost region (Xs) and region scattered into (IXs) */ 186 if (xs-sDist > 0) { 187 Xs = xs - sDist; 188 IXs = xs - sDist; 189 } else { 190 if (bx) Xs = xs - sDist; 191 else Xs = 0; 192 IXs = 0; 193 } 194 if (xe+sDist <= M) { 195 Xe = xe + sDist; 196 IXe = xe + sDist; 197 } else { 198 if (bx) Xe = xe + sDist; 199 else Xe = M; 200 IXe = M; 201 } 202 203 if (bx == DM_BOUNDARY_PERIODIC || bx == DM_BOUNDARY_MIRROR) { 204 Xs = xs - sDist; 205 Xe = xe + sDist; 206 IXs = xs - sDist; 207 IXe = xe + sDist; 208 } 209 210 /* allocate the base parallel and sequential vectors */ 211 dd->Nlocal = dof*x; 212 ierr = VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,NULL,&global);CHKERRQ(ierr); 213 dd->nlocal = dof*(Xe-Xs); 214 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,NULL,&local);CHKERRQ(ierr); 215 216 ierr = VecGetOwnershipRange(global,&start,NULL);CHKERRQ(ierr); 217 218 /* Create Global to Local Vector Scatter Context */ 219 /* global to local must retrieve ghost points */ 220 ierr = ISCreateStride(comm,dof*(IXe-IXs),dof*(IXs-Xs),1,&to);CHKERRQ(ierr); 221 222 ierr = PetscMalloc1(x+2*sDist,&idx);CHKERRQ(ierr); 223 ierr = PetscLogObjectMemory((PetscObject)da,(x+2*(sDist))*sizeof(PetscInt));CHKERRQ(ierr); 224 225 for (i=0; i<IXs-Xs; i++) idx[i] = -1; /* prepend with -1s if needed for ghosted case*/ 226 227 nn = IXs-Xs; 228 if (bx == DM_BOUNDARY_PERIODIC) { /* Handle all cases with periodic first */ 229 for (i=0; i<sDist; i++) { /* Left ghost points */ 230 if ((xs-sDist+i)>=0) idx[nn++] = xs-sDist+i; 231 else idx[nn++] = M+(xs-sDist+i); 232 } 233 234 for (i=0; i<x; i++) idx [nn++] = xs + i; /* Non-ghost points */ 235 236 for (i=0; i<sDist; i++) { /* Right ghost points */ 237 if ((xe+i)<M) idx [nn++] = xe+i; 238 else idx [nn++] = (xe+i) - M; 239 } 240 } else if (bx == DM_BOUNDARY_MIRROR) { /* Handle all cases with periodic first */ 241 for (i=0; i<(sDist); i++) { /* Left ghost points */ 242 if ((xs-sDist+i)>=0) idx[nn++] = xs-sDist+i; 243 else idx[nn++] = sDist - i; 244 } 245 246 for (i=0; i<x; i++) idx [nn++] = xs + i; /* Non-ghost points */ 247 248 for (i=0; i<(sDist); i++) { /* Right ghost points */ 249 if ((xe+i)<M) idx[nn++] = xe+i; 250 else idx[nn++] = M - (i + 1); 251 } 252 } else { /* Now do all cases with no periodicity */ 253 if (0 <= xs-sDist) { 254 for (i=0; i<sDist; i++) idx[nn++] = xs - sDist + i; 255 } else { 256 for (i=0; i<xs; i++) idx[nn++] = i; 257 } 258 259 for (i=0; i<x; i++) idx [nn++] = xs + i; 260 261 if ((xe+sDist)<=M) { 262 for (i=0; i<sDist; i++) idx[nn++]=xe+i; 263 } else { 264 for (i=xe; i<M; i++) idx[nn++]=i; 265 } 266 } 267 268 ierr = ISCreateBlock(comm,dof,nn-IXs+Xs,&idx[IXs-Xs],PETSC_USE_POINTER,&from);CHKERRQ(ierr); 269 ierr = VecScatterCreate(global,from,local,to,>ol);CHKERRQ(ierr); 270 ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)gtol);CHKERRQ(ierr); 271 ierr = ISDestroy(&to);CHKERRQ(ierr); 272 ierr = ISDestroy(&from);CHKERRQ(ierr); 273 ierr = VecDestroy(&local);CHKERRQ(ierr); 274 ierr = VecDestroy(&global);CHKERRQ(ierr); 275 276 dd->xs = dof*xs; dd->xe = dof*xe; dd->ys = 0; dd->ye = 1; dd->zs = 0; dd->ze = 1; 277 dd->Xs = dof*Xs; dd->Xe = dof*Xe; dd->Ys = 0; dd->Ye = 1; dd->Zs = 0; dd->Ze = 1; 278 279 dd->gtol = gtol; 280 dd->base = dof*xs; 281 da->ops->view = DMView_DA_1d; 282 283 /* 284 Set the local to global ordering in the global vector, this allows use 285 of VecSetValuesLocal(). 286 */ 287 for (i=0; i<Xe-IXe; i++) idx[nn++] = -1; /* pad with -1s if needed for ghosted case*/ 288 289 ierr = ISLocalToGlobalMappingCreate(comm,dof,nn,idx,PETSC_OWN_POINTER,&da->ltogmap);CHKERRQ(ierr); 290 ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)da->ltogmap);CHKERRQ(ierr); 291 292 PetscFunctionReturn(0); 293 } 294 295 296 #undef __FUNCT__ 297 #define __FUNCT__ "DMDACreate1d" 298 /*@C 299 DMDACreate1d - Creates an object that will manage the communication of one-dimensional 300 regular array data that is distributed across some processors. 301 302 Collective on MPI_Comm 303 304 Input Parameters: 305 + comm - MPI communicator 306 . bx - type of ghost cells at the boundary the array should have, if any. Use 307 DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, or DM_BOUNDARY_PERIODIC. 308 . M - global dimension of the array (use -M to indicate that it may be set to a different value 309 from the command line with -da_grid_x <M>) 310 . dof - number of degrees of freedom per node 311 . s - stencil width 312 - lx - array containing number of nodes in the X direction on each processor, 313 or NULL. If non-null, must be of length as the number of processes in the MPI_Comm. 314 315 Output Parameter: 316 . da - the resulting distributed array object 317 318 Options Database Key: 319 + -dm_view - Calls DMView() at the conclusion of DMDACreate1d() 320 . -da_grid_x <nx> - number of grid points in x direction; can set if M < 0 321 . -da_refine_x <rx> - refinement factor 322 - -da_refine <n> - refine the DMDA n times before creating it, if M < 0 323 324 Level: beginner 325 326 Notes: 327 The array data itself is NOT stored in the DMDA, it is stored in Vec objects; 328 The appropriate vector objects can be obtained with calls to DMCreateGlobalVector() 329 and DMCreateLocalVector() and calls to VecDuplicate() if more are needed. 330 331 .keywords: distributed array, create, one-dimensional 332 333 .seealso: DMDestroy(), DMView(), DMDACreate2d(), DMDACreate3d(), DMGlobalToLocalBegin(), DMDASetRefinementFactor(), 334 DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMLocalToLocalBegin(), DMLocalToLocalEnd(), DMDAGetRefinementFactor(), 335 DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges() 336 337 @*/ 338 PetscErrorCode DMDACreate1d(MPI_Comm comm, DMBoundaryType bx, PetscInt M, PetscInt dof, PetscInt s, const PetscInt lx[], DM *da) 339 { 340 PetscErrorCode ierr; 341 PetscMPIInt size; 342 343 PetscFunctionBegin; 344 ierr = DMDACreate(comm, da);CHKERRQ(ierr); 345 ierr = DMSetDimension(*da, 1);CHKERRQ(ierr); 346 ierr = DMDASetSizes(*da, M, 1, 1);CHKERRQ(ierr); 347 ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 348 ierr = DMDASetNumProcs(*da, size, PETSC_DECIDE, PETSC_DECIDE);CHKERRQ(ierr); 349 ierr = DMDASetBoundaryType(*da, bx, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE);CHKERRQ(ierr); 350 ierr = DMDASetDof(*da, dof);CHKERRQ(ierr); 351 ierr = DMDASetStencilWidth(*da, s);CHKERRQ(ierr); 352 ierr = DMDASetOwnershipRanges(*da, lx, NULL, NULL);CHKERRQ(ierr); 353 /* This violates the behavior for other classes, but right now users expect negative dimensions to be handled this way */ 354 ierr = DMSetFromOptions(*da);CHKERRQ(ierr); 355 ierr = DMSetUp(*da);CHKERRQ(ierr); 356 PetscFunctionReturn(0); 357 } 358