Actual source code: da2.c
petsc-3.8.3 2017-12-09
2: #include <petsc/private/dmdaimpl.h>
3: #include <petscdraw.h>
5: static PetscErrorCode DMView_DA_2d(DM da,PetscViewer viewer)
6: {
8: PetscMPIInt rank;
9: PetscBool iascii,isdraw,isglvis,isbinary;
10: DM_DA *dd = (DM_DA*)da->data;
11: #if defined(PETSC_HAVE_MATLAB_ENGINE)
12: PetscBool ismatlab;
13: #endif
16: MPI_Comm_rank(PetscObjectComm((PetscObject)da),&rank);
18: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
19: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
20: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERGLVIS,&isglvis);
21: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
22: #if defined(PETSC_HAVE_MATLAB_ENGINE)
23: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERMATLAB,&ismatlab);
24: #endif
25: if (iascii) {
26: PetscViewerFormat format;
28: PetscViewerGetFormat(viewer, &format);
29: if (format != PETSC_VIEWER_ASCII_VTK && format != PETSC_VIEWER_ASCII_VTK_CELL && format != PETSC_VIEWER_ASCII_GLVIS) {
30: DMDALocalInfo info;
31: DMDAGetLocalInfo(da,&info);
32: PetscViewerASCIIPushSynchronized(viewer);
33: 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);
34: 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);
35: PetscViewerFlush(viewer);
36: PetscViewerASCIIPopSynchronized(viewer);
37: } else if (format == PETSC_VIEWER_ASCII_GLVIS) {
38: DMView_DA_GLVis(da,viewer);
39: } else {
40: DMView_DA_VTK(da,viewer);
41: }
42: } else if (isdraw) {
43: PetscDraw draw;
44: double ymin = -1*dd->s-1,ymax = dd->N+dd->s;
45: double xmin = -1*dd->s-1,xmax = dd->M+dd->s;
46: double x,y;
47: PetscInt base;
48: const PetscInt *idx;
49: char node[10];
50: PetscBool isnull;
52: PetscViewerDrawGetDraw(viewer,0,&draw);
53: PetscDrawIsNull(draw,&isnull);
54: if (isnull) return(0);
56: PetscDrawCheckResizedWindow(draw);
57: PetscDrawClear(draw);
58: PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);
60: PetscDrawCollectiveBegin(draw);
61: /* first processor draw all node lines */
62: if (!rank) {
63: ymin = 0.0; ymax = dd->N - 1;
64: for (xmin=0; xmin<dd->M; xmin++) {
65: PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK);
66: }
67: xmin = 0.0; xmax = dd->M - 1;
68: for (ymin=0; ymin<dd->N; ymin++) {
69: PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);
70: }
71: }
72: PetscDrawCollectiveEnd(draw);
73: PetscDrawFlush(draw);
74: PetscDrawPause(draw);
76: PetscDrawCollectiveBegin(draw);
77: /* draw my box */
78: xmin = dd->xs/dd->w; xmax =(dd->xe-1)/dd->w; ymin = dd->ys; ymax = dd->ye - 1;
79: PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);
80: PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);
81: PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);
82: PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);
83: /* put in numbers */
84: base = (dd->base)/dd->w;
85: for (y=ymin; y<=ymax; y++) {
86: for (x=xmin; x<=xmax; x++) {
87: PetscSNPrintf(node,sizeof(node),"%d",(int)base++);
88: PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node);
89: }
90: }
91: PetscDrawCollectiveEnd(draw);
92: PetscDrawFlush(draw);
93: PetscDrawPause(draw);
95: PetscDrawCollectiveBegin(draw);
96: /* overlay ghost numbers, useful for error checking */
97: ISLocalToGlobalMappingGetBlockIndices(da->ltogmap,&idx);
98: base = 0; xmin = dd->Xs; xmax = dd->Xe; ymin = dd->Ys; ymax = dd->Ye;
99: for (y=ymin; y<ymax; y++) {
100: for (x=xmin; x<xmax; x++) {
101: if ((base % dd->w) == 0) {
102: PetscSNPrintf(node,sizeof(node),"%d",(int)(idx[base/dd->w]));
103: PetscDrawString(draw,x/dd->w,y,PETSC_DRAW_BLUE,node);
104: }
105: base++;
106: }
107: }
108: ISLocalToGlobalMappingRestoreBlockIndices(da->ltogmap,&idx);
109: PetscDrawCollectiveEnd(draw);
110: PetscDrawFlush(draw);
111: PetscDrawPause(draw);
112: PetscDrawSave(draw);
113: } else if (isglvis) {
114: DMView_DA_GLVis(da,viewer);
115: } else if (isbinary) {
116: DMView_DA_Binary(da,viewer);
117: #if defined(PETSC_HAVE_MATLAB_ENGINE)
118: } else if (ismatlab) {
119: DMView_DA_Matlab(da,viewer);
120: #endif
121: }
122: return(0);
123: }
125: /*
126: M is number of grid points
127: m is number of processors
129: */
130: PetscErrorCode DMDASplitComm2d(MPI_Comm comm,PetscInt M,PetscInt N,PetscInt sw,MPI_Comm *outcomm)
131: {
133: PetscInt m,n = 0,x = 0,y = 0;
134: PetscMPIInt size,csize,rank;
137: MPI_Comm_size(comm,&size);
138: MPI_Comm_rank(comm,&rank);
140: csize = 4*size;
141: do {
142: if (csize % 4) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cannot split communicator of size %d tried %d %D %D",size,csize,x,y);
143: csize = csize/4;
145: m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)csize)/((PetscReal)N)));
146: if (!m) m = 1;
147: while (m > 0) {
148: n = csize/m;
149: if (m*n == csize) break;
150: m--;
151: }
152: if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
154: x = M/m + ((M % m) > ((csize-1) % m));
155: y = (N + (csize-1)/m)/n;
156: } while ((x < 4 || y < 4) && csize > 1);
157: if (size != csize) {
158: MPI_Group entire_group,sub_group;
159: PetscMPIInt i,*groupies;
161: MPI_Comm_group(comm,&entire_group);
162: PetscMalloc1(csize,&groupies);
163: for (i=0; i<csize; i++) {
164: groupies[i] = (rank/csize)*csize + i;
165: }
166: MPI_Group_incl(entire_group,csize,groupies,&sub_group);
167: PetscFree(groupies);
168: MPI_Comm_create(comm,sub_group,outcomm);
169: MPI_Group_free(&entire_group);
170: MPI_Group_free(&sub_group);
171: PetscInfo1(0,"DMDASplitComm2d:Creating redundant coarse problems of size %d\n",csize);
172: } else {
173: *outcomm = comm;
174: }
175: return(0);
176: }
178: #if defined(new)
179: /*
180: DMDAGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix where local
181: function lives on a DMDA
183: y ~= (F(u + ha) - F(u))/h,
184: where F = nonlinear function, as set by SNESSetFunction()
185: u = current iterate
186: h = difference interval
187: */
188: PetscErrorCode DMDAGetDiagonal_MFFD(DM da,Vec U,Vec a)
189: {
190: PetscScalar h,*aa,*ww,v;
191: PetscReal epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
193: PetscInt gI,nI;
194: MatStencil stencil;
195: DMDALocalInfo info;
198: (*ctx->func)(0,U,a,ctx->funcctx);
199: (*ctx->funcisetbase)(U,ctx->funcctx);
201: VecGetArray(U,&ww);
202: VecGetArray(a,&aa);
204: nI = 0;
205: h = ww[gI];
206: if (h == 0.0) h = 1.0;
207: if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin;
208: else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
209: h *= epsilon;
211: ww[gI] += h;
212: (*ctx->funci)(i,w,&v,ctx->funcctx);
213: aa[nI] = (v - aa[nI])/h;
214: ww[gI] -= h;
215: nI++;
217: VecRestoreArray(U,&ww);
218: VecRestoreArray(a,&aa);
219: return(0);
220: }
221: #endif
223: PetscErrorCode DMSetUp_DA_2D(DM da)
224: {
225: DM_DA *dd = (DM_DA*)da->data;
226: const PetscInt M = dd->M;
227: const PetscInt N = dd->N;
228: PetscInt m = dd->m;
229: PetscInt n = dd->n;
230: const PetscInt dof = dd->w;
231: const PetscInt s = dd->s;
232: DMBoundaryType bx = dd->bx;
233: DMBoundaryType by = dd->by;
234: DMDAStencilType stencil_type = dd->stencil_type;
235: PetscInt *lx = dd->lx;
236: PetscInt *ly = dd->ly;
237: MPI_Comm comm;
238: PetscMPIInt rank,size;
239: PetscInt xs,xe,ys,ye,x,y,Xs,Xe,Ys,Ye,IXs,IXe,IYs,IYe;
240: PetscInt up,down,left,right,i,n0,n1,n2,n3,n5,n6,n7,n8,*idx,nn;
241: PetscInt xbase,*bases,*ldims,j,x_t,y_t,s_t,base,count;
242: PetscInt s_x,s_y; /* s proportionalized to w */
243: PetscInt sn0 = 0,sn2 = 0,sn6 = 0,sn8 = 0;
244: Vec local,global;
245: VecScatter gtol;
246: IS to,from;
247: PetscErrorCode ierr;
250: if (stencil_type == DMDA_STENCIL_BOX && (bx == DM_BOUNDARY_MIRROR || by == DM_BOUNDARY_MIRROR)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Mirror boundary and box stencil");
251: PetscObjectGetComm((PetscObject)da,&comm);
252: #if !defined(PETSC_USE_64BIT_INDICES)
253: if (((PetscInt64) M)*((PetscInt64) N)*((PetscInt64) dof) > (PetscInt64) PETSC_MPI_INT_MAX) SETERRQ3(comm,PETSC_ERR_INT_OVERFLOW,"Mesh of %D by %D by %D (dof) is too large for 32 bit indices",M,N,dof);
254: #endif
256: if (dof < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have 1 or more degrees of freedom per node: %D",dof);
257: if (s < 0) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Stencil width cannot be negative: %D",s);
259: MPI_Comm_size(comm,&size);
260: MPI_Comm_rank(comm,&rank);
262: dd->p = 1;
263: if (m != PETSC_DECIDE) {
264: if (m < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m);
265: else if (m > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size);
266: }
267: if (n != PETSC_DECIDE) {
268: if (n < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n);
269: else if (n > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size);
270: }
272: if (m == PETSC_DECIDE || n == PETSC_DECIDE) {
273: if (n != PETSC_DECIDE) {
274: m = size/n;
275: } else if (m != PETSC_DECIDE) {
276: n = size/m;
277: } else {
278: /* try for squarish distribution */
279: m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)N)));
280: if (!m) m = 1;
281: while (m > 0) {
282: n = size/m;
283: if (m*n == size) break;
284: m--;
285: }
286: if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
287: }
288: if (m*n != size) SETERRQ(comm,PETSC_ERR_PLIB,"Unable to create partition, check the size of the communicator and input m and n ");
289: } else if (m*n != size) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition");
291: if (M < m) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m);
292: if (N < n) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n);
294: /*
295: Determine locally owned region
296: xs is the first local node number, x is the number of local nodes
297: */
298: if (!lx) {
299: PetscMalloc1(m, &dd->lx);
300: lx = dd->lx;
301: for (i=0; i<m; i++) {
302: lx[i] = M/m + ((M % m) > i);
303: }
304: }
305: x = lx[rank % m];
306: xs = 0;
307: for (i=0; i<(rank % m); i++) {
308: xs += lx[i];
309: }
310: #if defined(PETSC_USE_DEBUG)
311: left = xs;
312: for (i=(rank % m); i<m; i++) {
313: left += lx[i];
314: }
315: if (left != M) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Sum of lx across processors not equal to M: %D %D",left,M);
316: #endif
318: /*
319: Determine locally owned region
320: ys is the first local node number, y is the number of local nodes
321: */
322: if (!ly) {
323: PetscMalloc1(n, &dd->ly);
324: ly = dd->ly;
325: for (i=0; i<n; i++) {
326: ly[i] = N/n + ((N % n) > i);
327: }
328: }
329: y = ly[rank/m];
330: ys = 0;
331: for (i=0; i<(rank/m); i++) {
332: ys += ly[i];
333: }
334: #if defined(PETSC_USE_DEBUG)
335: left = ys;
336: for (i=(rank/m); i<n; i++) {
337: left += ly[i];
338: }
339: if (left != N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Sum of ly across processors not equal to N: %D %D",left,N);
340: #endif
342: /*
343: check if the scatter requires more than one process neighbor or wraps around
344: the domain more than once
345: */
346: 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);
347: if ((y < s) && ((n > 1) || (by == DM_BOUNDARY_PERIODIC))) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local y-width of domain y %D is smaller than stencil width s %D",y,s);
348: xe = xs + x;
349: ye = ys + y;
351: /* determine ghost region (Xs) and region scattered into (IXs) */
352: if (xs-s > 0) {
353: Xs = xs - s; IXs = xs - s;
354: } else {
355: if (bx) {
356: Xs = xs - s;
357: } else {
358: Xs = 0;
359: }
360: IXs = 0;
361: }
362: if (xe+s <= M) {
363: Xe = xe + s; IXe = xe + s;
364: } else {
365: if (bx) {
366: Xs = xs - s; Xe = xe + s;
367: } else {
368: Xe = M;
369: }
370: IXe = M;
371: }
373: if (bx == DM_BOUNDARY_PERIODIC || bx == DM_BOUNDARY_MIRROR) {
374: IXs = xs - s;
375: IXe = xe + s;
376: Xs = xs - s;
377: Xe = xe + s;
378: }
380: if (ys-s > 0) {
381: Ys = ys - s; IYs = ys - s;
382: } else {
383: if (by) {
384: Ys = ys - s;
385: } else {
386: Ys = 0;
387: }
388: IYs = 0;
389: }
390: if (ye+s <= N) {
391: Ye = ye + s; IYe = ye + s;
392: } else {
393: if (by) {
394: Ye = ye + s;
395: } else {
396: Ye = N;
397: }
398: IYe = N;
399: }
401: if (by == DM_BOUNDARY_PERIODIC || by == DM_BOUNDARY_MIRROR) {
402: IYs = ys - s;
403: IYe = ye + s;
404: Ys = ys - s;
405: Ye = ye + s;
406: }
408: /* stencil length in each direction */
409: s_x = s;
410: s_y = s;
412: /* determine starting point of each processor */
413: nn = x*y;
414: PetscMalloc2(size+1,&bases,size,&ldims);
415: MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);
416: bases[0] = 0;
417: for (i=1; i<=size; i++) {
418: bases[i] = ldims[i-1];
419: }
420: for (i=1; i<=size; i++) {
421: bases[i] += bases[i-1];
422: }
423: base = bases[rank]*dof;
425: /* allocate the base parallel and sequential vectors */
426: dd->Nlocal = x*y*dof;
427: VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,NULL,&global);
428: dd->nlocal = (Xe-Xs)*(Ye-Ys)*dof;
429: VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,NULL,&local);
431: /* generate appropriate vector scatters */
432: /* local to global inserts non-ghost point region into global */
433: PetscMalloc1((IXe-IXs)*(IYe-IYs),&idx);
434: left = xs - Xs; right = left + x;
435: down = ys - Ys; up = down + y;
436: count = 0;
437: for (i=down; i<up; i++) {
438: for (j=left; j<right; j++) {
439: idx[count++] = i*(Xe-Xs) + j;
440: }
441: }
443: /* global to local must include ghost points within the domain,
444: but not ghost points outside the domain that aren't periodic */
445: if (stencil_type == DMDA_STENCIL_BOX) {
446: left = IXs - Xs; right = left + (IXe-IXs);
447: down = IYs - Ys; up = down + (IYe-IYs);
448: count = 0;
449: for (i=down; i<up; i++) {
450: for (j=left; j<right; j++) {
451: idx[count++] = j + i*(Xe-Xs);
452: }
453: }
454: ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);
456: } else {
457: /* must drop into cross shape region */
458: /* ---------|
459: | top |
460: |--- ---| up
461: | middle |
462: | |
463: ---- ---- down
464: | bottom |
465: -----------
466: Xs xs xe Xe */
467: left = xs - Xs; right = left + x;
468: down = ys - Ys; up = down + y;
469: count = 0;
470: /* bottom */
471: for (i=(IYs-Ys); i<down; i++) {
472: for (j=left; j<right; j++) {
473: idx[count++] = j + i*(Xe-Xs);
474: }
475: }
476: /* middle */
477: for (i=down; i<up; i++) {
478: for (j=(IXs-Xs); j<(IXe-Xs); j++) {
479: idx[count++] = j + i*(Xe-Xs);
480: }
481: }
482: /* top */
483: for (i=up; i<up+IYe-ye; i++) {
484: for (j=left; j<right; j++) {
485: idx[count++] = j + i*(Xe-Xs);
486: }
487: }
488: ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);
489: }
492: /* determine who lies on each side of us stored in n6 n7 n8
493: n3 n5
494: n0 n1 n2
495: */
497: /* Assume the Non-Periodic Case */
498: n1 = rank - m;
499: if (rank % m) {
500: n0 = n1 - 1;
501: } else {
502: n0 = -1;
503: }
504: if ((rank+1) % m) {
505: n2 = n1 + 1;
506: n5 = rank + 1;
507: n8 = rank + m + 1; if (n8 >= m*n) n8 = -1;
508: } else {
509: n2 = -1; n5 = -1; n8 = -1;
510: }
511: if (rank % m) {
512: n3 = rank - 1;
513: n6 = n3 + m; if (n6 >= m*n) n6 = -1;
514: } else {
515: n3 = -1; n6 = -1;
516: }
517: n7 = rank + m; if (n7 >= m*n) n7 = -1;
519: if (bx == DM_BOUNDARY_PERIODIC && by == DM_BOUNDARY_PERIODIC) {
520: /* Modify for Periodic Cases */
521: /* Handle all four corners */
522: if ((n6 < 0) && (n7 < 0) && (n3 < 0)) n6 = m-1;
523: if ((n8 < 0) && (n7 < 0) && (n5 < 0)) n8 = 0;
524: if ((n2 < 0) && (n5 < 0) && (n1 < 0)) n2 = size-m;
525: if ((n0 < 0) && (n3 < 0) && (n1 < 0)) n0 = size-1;
527: /* Handle Top and Bottom Sides */
528: if (n1 < 0) n1 = rank + m * (n-1);
529: if (n7 < 0) n7 = rank - m * (n-1);
530: if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
531: if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
532: if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
533: if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
535: /* Handle Left and Right Sides */
536: if (n3 < 0) n3 = rank + (m-1);
537: if (n5 < 0) n5 = rank - (m-1);
538: if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
539: if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
540: if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
541: if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
542: } else if (by == DM_BOUNDARY_PERIODIC) { /* Handle Top and Bottom Sides */
543: if (n1 < 0) n1 = rank + m * (n-1);
544: if (n7 < 0) n7 = rank - m * (n-1);
545: if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
546: if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
547: if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
548: if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
549: } else if (bx == DM_BOUNDARY_PERIODIC) { /* Handle Left and Right Sides */
550: if (n3 < 0) n3 = rank + (m-1);
551: if (n5 < 0) n5 = rank - (m-1);
552: if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
553: if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
554: if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
555: if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
556: }
558: PetscMalloc1(9,&dd->neighbors);
560: dd->neighbors[0] = n0;
561: dd->neighbors[1] = n1;
562: dd->neighbors[2] = n2;
563: dd->neighbors[3] = n3;
564: dd->neighbors[4] = rank;
565: dd->neighbors[5] = n5;
566: dd->neighbors[6] = n6;
567: dd->neighbors[7] = n7;
568: dd->neighbors[8] = n8;
570: if (stencil_type == DMDA_STENCIL_STAR) {
571: /* save corner processor numbers */
572: sn0 = n0; sn2 = n2; sn6 = n6; sn8 = n8;
573: n0 = n2 = n6 = n8 = -1;
574: }
576: PetscMalloc1((Xe-Xs)*(Ye-Ys),&idx);
578: nn = 0;
579: xbase = bases[rank];
580: for (i=1; i<=s_y; i++) {
581: if (n0 >= 0) { /* left below */
582: x_t = lx[n0 % m];
583: y_t = ly[(n0/m)];
584: s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
585: for (j=0; j<s_x; j++) idx[nn++] = s_t++;
586: }
588: if (n1 >= 0) { /* directly below */
589: x_t = x;
590: y_t = ly[(n1/m)];
591: s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
592: for (j=0; j<x_t; j++) idx[nn++] = s_t++;
593: } else if (by == DM_BOUNDARY_MIRROR) {
594: for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(s_y - i + 1) + j;
595: }
597: if (n2 >= 0) { /* right below */
598: x_t = lx[n2 % m];
599: y_t = ly[(n2/m)];
600: s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
601: for (j=0; j<s_x; j++) idx[nn++] = s_t++;
602: }
603: }
605: for (i=0; i<y; i++) {
606: if (n3 >= 0) { /* directly left */
607: x_t = lx[n3 % m];
608: /* y_t = y; */
609: s_t = bases[n3] + (i+1)*x_t - s_x;
610: for (j=0; j<s_x; j++) idx[nn++] = s_t++;
611: } else if (bx == DM_BOUNDARY_MIRROR) {
612: for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*i + s_x - j;
613: }
615: for (j=0; j<x; j++) idx[nn++] = xbase++; /* interior */
617: if (n5 >= 0) { /* directly right */
618: x_t = lx[n5 % m];
619: /* y_t = y; */
620: s_t = bases[n5] + (i)*x_t;
621: for (j=0; j<s_x; j++) idx[nn++] = s_t++;
622: } else if (bx == DM_BOUNDARY_MIRROR) {
623: for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*(i + 1) - 2 - j;
624: }
625: }
627: for (i=1; i<=s_y; i++) {
628: if (n6 >= 0) { /* left above */
629: x_t = lx[n6 % m];
630: /* y_t = ly[(n6/m)]; */
631: s_t = bases[n6] + (i)*x_t - s_x;
632: for (j=0; j<s_x; j++) idx[nn++] = s_t++;
633: }
635: if (n7 >= 0) { /* directly above */
636: x_t = x;
637: /* y_t = ly[(n7/m)]; */
638: s_t = bases[n7] + (i-1)*x_t;
639: for (j=0; j<x_t; j++) idx[nn++] = s_t++;
640: } else if (by == DM_BOUNDARY_MIRROR) {
641: for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(y - i - 1) + j;
642: }
644: if (n8 >= 0) { /* right above */
645: x_t = lx[n8 % m];
646: /* y_t = ly[(n8/m)]; */
647: s_t = bases[n8] + (i-1)*x_t;
648: for (j=0; j<s_x; j++) idx[nn++] = s_t++;
649: }
650: }
652: ISCreateBlock(comm,dof,nn,idx,PETSC_USE_POINTER,&from);
653: VecScatterCreate(global,from,local,to,>ol);
654: PetscLogObjectParent((PetscObject)da,(PetscObject)gtol);
655: ISDestroy(&to);
656: ISDestroy(&from);
658: if (stencil_type == DMDA_STENCIL_STAR) {
659: n0 = sn0; n2 = sn2; n6 = sn6; n8 = sn8;
660: }
662: if (((stencil_type == DMDA_STENCIL_STAR) || (bx && bx != DM_BOUNDARY_PERIODIC) || (by && by != DM_BOUNDARY_PERIODIC))) {
663: /*
664: Recompute the local to global mappings, this time keeping the
665: information about the cross corner processor numbers and any ghosted
666: but not periodic indices.
667: */
668: nn = 0;
669: xbase = bases[rank];
670: for (i=1; i<=s_y; i++) {
671: if (n0 >= 0) { /* left below */
672: x_t = lx[n0 % m];
673: y_t = ly[(n0/m)];
674: s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
675: for (j=0; j<s_x; j++) idx[nn++] = s_t++;
676: } else if (xs-Xs > 0 && ys-Ys > 0) {
677: for (j=0; j<s_x; j++) idx[nn++] = -1;
678: }
679: if (n1 >= 0) { /* directly below */
680: x_t = x;
681: y_t = ly[(n1/m)];
682: s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
683: for (j=0; j<x_t; j++) idx[nn++] = s_t++;
684: } else if (ys-Ys > 0) {
685: if (by == DM_BOUNDARY_MIRROR) {
686: for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(s_y - i + 1) + j;
687: } else {
688: for (j=0; j<x; j++) idx[nn++] = -1;
689: }
690: }
691: if (n2 >= 0) { /* right below */
692: x_t = lx[n2 % m];
693: y_t = ly[(n2/m)];
694: s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
695: for (j=0; j<s_x; j++) idx[nn++] = s_t++;
696: } else if (Xe-xe> 0 && ys-Ys > 0) {
697: for (j=0; j<s_x; j++) idx[nn++] = -1;
698: }
699: }
701: for (i=0; i<y; i++) {
702: if (n3 >= 0) { /* directly left */
703: x_t = lx[n3 % m];
704: /* y_t = y; */
705: s_t = bases[n3] + (i+1)*x_t - s_x;
706: for (j=0; j<s_x; j++) idx[nn++] = s_t++;
707: } else if (xs-Xs > 0) {
708: if (bx == DM_BOUNDARY_MIRROR) {
709: for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*i + s_x - j;
710: } else {
711: for (j=0; j<s_x; j++) idx[nn++] = -1;
712: }
713: }
715: for (j=0; j<x; j++) idx[nn++] = xbase++; /* interior */
717: if (n5 >= 0) { /* directly right */
718: x_t = lx[n5 % m];
719: /* y_t = y; */
720: s_t = bases[n5] + (i)*x_t;
721: for (j=0; j<s_x; j++) idx[nn++] = s_t++;
722: } else if (Xe-xe > 0) {
723: if (bx == DM_BOUNDARY_MIRROR) {
724: for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*(i + 1) - 2 - j;
725: } else {
726: for (j=0; j<s_x; j++) idx[nn++] = -1;
727: }
728: }
729: }
731: for (i=1; i<=s_y; i++) {
732: if (n6 >= 0) { /* left above */
733: x_t = lx[n6 % m];
734: /* y_t = ly[(n6/m)]; */
735: s_t = bases[n6] + (i)*x_t - s_x;
736: for (j=0; j<s_x; j++) idx[nn++] = s_t++;
737: } else if (xs-Xs > 0 && Ye-ye > 0) {
738: for (j=0; j<s_x; j++) idx[nn++] = -1;
739: }
740: if (n7 >= 0) { /* directly above */
741: x_t = x;
742: /* y_t = ly[(n7/m)]; */
743: s_t = bases[n7] + (i-1)*x_t;
744: for (j=0; j<x_t; j++) idx[nn++] = s_t++;
745: } else if (Ye-ye > 0) {
746: if (by == DM_BOUNDARY_MIRROR) {
747: for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(y - i - 1) + j;
748: } else {
749: for (j=0; j<x; j++) idx[nn++] = -1;
750: }
751: }
752: if (n8 >= 0) { /* right above */
753: x_t = lx[n8 % m];
754: /* y_t = ly[(n8/m)]; */
755: s_t = bases[n8] + (i-1)*x_t;
756: for (j=0; j<s_x; j++) idx[nn++] = s_t++;
757: } else if (Xe-xe > 0 && Ye-ye > 0) {
758: for (j=0; j<s_x; j++) idx[nn++] = -1;
759: }
760: }
761: }
762: /*
763: Set the local to global ordering in the global vector, this allows use
764: of VecSetValuesLocal().
765: */
766: ISLocalToGlobalMappingCreate(comm,dof,nn,idx,PETSC_OWN_POINTER,&da->ltogmap);
767: PetscLogObjectParent((PetscObject)da,(PetscObject)da->ltogmap);
769: PetscFree2(bases,ldims);
770: dd->m = m; dd->n = n;
771: /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */
772: dd->xs = xs*dof; dd->xe = xe*dof; dd->ys = ys; dd->ye = ye; dd->zs = 0; dd->ze = 1;
773: dd->Xs = Xs*dof; dd->Xe = Xe*dof; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = 0; dd->Ze = 1;
775: VecDestroy(&local);
776: VecDestroy(&global);
778: dd->gtol = gtol;
779: dd->base = base;
780: da->ops->view = DMView_DA_2d;
781: dd->ltol = NULL;
782: dd->ao = NULL;
783: return(0);
784: }
786: /*@C
787: DMDACreate2d - Creates an object that will manage the communication of two-dimensional
788: regular array data that is distributed across some processors.
790: Collective on MPI_Comm
792: Input Parameters:
793: + comm - MPI communicator
794: . bx,by - type of ghost nodes the array have.
795: Use one of DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_PERIODIC.
796: . stencil_type - stencil type. Use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
797: . M,N - global dimension in each direction of the array
798: . m,n - corresponding number of processors in each dimension
799: (or PETSC_DECIDE to have calculated)
800: . dof - number of degrees of freedom per node
801: . s - stencil width
802: - lx, ly - arrays containing the number of nodes in each cell along
803: the x and y coordinates, or NULL. If non-null, these
804: must be of length as m and n, and the corresponding
805: m and n cannot be PETSC_DECIDE. The sum of the lx[] entries
806: must be M, and the sum of the ly[] entries must be N.
808: Output Parameter:
809: . da - the resulting distributed array object
811: Options Database Key:
812: + -dm_view - Calls DMView() at the conclusion of DMDACreate2d()
813: . -da_grid_x <nx> - number of grid points in x direction, if M < 0
814: . -da_grid_y <ny> - number of grid points in y direction, if N < 0
815: . -da_processors_x <nx> - number of processors in x direction
816: . -da_processors_y <ny> - number of processors in y direction
817: . -da_refine_x <rx> - refinement ratio in x direction
818: . -da_refine_y <ry> - refinement ratio in y direction
819: - -da_refine <n> - refine the DMDA n times before creating, if M or N < 0
822: Level: beginner
824: Notes:
825: The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the
826: standard 5-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes
827: the standard 9-pt stencil.
829: The array data itself is NOT stored in the DMDA, it is stored in Vec objects;
830: The appropriate vector objects can be obtained with calls to DMCreateGlobalVector()
831: and DMCreateLocalVector() and calls to VecDuplicate() if more are needed.
833: You must call DMSetUp() after this call before using this DM.
835: If you wish to use the options database to change values in the DMDA call DMSetFromOptions() after this call
836: but before DMSetUp().
838: .keywords: distributed array, create, two-dimensional
840: .seealso: DMDestroy(), DMView(), DMDACreate1d(), DMDACreate3d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(),
841: DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMLocalToLocalBegin(), DMLocalToLocalEnd(), DMDASetRefinementFactor(),
842: DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges()
844: @*/
846: PetscErrorCode DMDACreate2d(MPI_Comm comm,DMBoundaryType bx,DMBoundaryType by,DMDAStencilType stencil_type,
847: PetscInt M,PetscInt N,PetscInt m,PetscInt n,PetscInt dof,PetscInt s,const PetscInt lx[],const PetscInt ly[],DM *da)
848: {
852: DMDACreate(comm, da);
853: DMSetDimension(*da, 2);
854: DMDASetSizes(*da, M, N, 1);
855: DMDASetNumProcs(*da, m, n, PETSC_DECIDE);
856: DMDASetBoundaryType(*da, bx, by, DM_BOUNDARY_NONE);
857: DMDASetDof(*da, dof);
858: DMDASetStencilType(*da, stencil_type);
859: DMDASetStencilWidth(*da, s);
860: DMDASetOwnershipRanges(*da, lx, ly, NULL);
861: return(0);
862: }