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