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