xref: /petsc/src/dm/impls/da/gr1.c (revision 5edf65166e22df6feda2ed6bd8a6e348419dee17)
1 
2 /*
3    Plots vectors obtained with DMDACreate1d()
4 */
5 
6 #include <petscdmda.h>      /*I  "petscdmda.h"   I*/
7 
8 #undef __FUNCT__
9 #define __FUNCT__ "DMDASetUniformCoordinates"
10 /*@
11     DMDASetUniformCoordinates - Sets a DMDA coordinates to be a uniform grid
12 
13   Collective on DMDA
14 
15   Input Parameters:
16 +  da - the distributed array object
17 .  xmin,xmax - extremes in the x direction
18 .  ymin,ymax - extremes in the y direction (use PETSC_NULL for 1 dimensional problems)
19 -  zmin,zmax - extremes in the z direction (use PETSC_NULL for 1 or 2 dimensional problems)
20 
21   Level: beginner
22 
23 .seealso: DMSetCoordinates(), DMGetCoordinates(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
24 
25 @*/
26 PetscErrorCode  DMDASetUniformCoordinates(DM da,PetscReal xmin,PetscReal xmax,PetscReal ymin,PetscReal ymax,PetscReal zmin,PetscReal zmax)
27 {
28   MPI_Comm         comm;
29   PetscSection     section;
30   DM               cda;
31   DMDABoundaryType bx,by,bz;
32   Vec              xcoor;
33   PetscScalar      *coors;
34   PetscReal        hx,hy,hz_;
35   PetscInt         i,j,k,M,N,P,istart,isize,jstart,jsize,kstart,ksize,dim,cnt;
36   PetscErrorCode   ierr;
37 
38   PetscFunctionBegin;
39   PetscValidHeaderSpecific(da,DM_CLASSID,1);
40   ierr = DMDAGetInfo(da,&dim,&M,&N,&P,0,0,0,0,0,&bx,&by,&bz,0);CHKERRQ(ierr);
41   if (xmax <= xmin) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_INCOMP,"xmax must be larger than xmin %G %G",xmin,xmax);
42   if ((ymax <= ymin) && (dim > 1)) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_INCOMP,"ymax must be larger than ymin %G %G",ymin,ymax);
43   if ((zmax <= zmin) && (dim > 2)) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_INCOMP,"zmax must be larger than zmin %G %G",zmin,zmax);
44   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
45   ierr = DMGetDefaultSection(da,&section);CHKERRQ(ierr);
46   ierr = DMDAGetCorners(da,&istart,&jstart,&kstart,&isize,&jsize,&ksize);CHKERRQ(ierr);
47   ierr = DMGetCoordinateDM(da, &cda);CHKERRQ(ierr);
48   if (section) {
49     /* This would be better as a vector, but this is compatible */
50     PetscInt numComp[3]      = {1, 1, 1};
51     PetscInt numVertexDof[3] = {1, 1, 1};
52 
53     ierr = DMDASetFieldName(cda, 0, "x");CHKERRQ(ierr);
54     if (dim > 1) {ierr = DMDASetFieldName(cda, 1, "y");CHKERRQ(ierr);}
55     if (dim > 2) {ierr = DMDASetFieldName(cda, 2, "z");CHKERRQ(ierr);}
56     ierr = DMDACreateSection(cda, numComp, numVertexDof, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
57   }
58   ierr = DMCreateGlobalVector(cda, &xcoor);CHKERRQ(ierr);
59   if (section) {
60     PetscSection csection;
61     PetscInt     vStart, vEnd;
62 
63     ierr = DMGetDefaultGlobalSection(cda,&csection);CHKERRQ(ierr);
64     ierr = VecGetArray(xcoor,&coors);CHKERRQ(ierr);
65     ierr = DMDAGetHeightStratum(da, dim, &vStart, &vEnd);CHKERRQ(ierr);
66     if (bx == DMDA_BOUNDARY_PERIODIC) hx  = (xmax-xmin)/(M+1);
67     else                              hx  = (xmax-xmin)/(M ? M : 1);
68     if (by == DMDA_BOUNDARY_PERIODIC) hy  = (ymax-ymin)/(N+1);
69     else                              hy  = (ymax-ymin)/(N ? N : 1);
70     if (bz == DMDA_BOUNDARY_PERIODIC) hz_ = (zmax-zmin)/(P+1);
71     else                              hz_ = (zmax-zmin)/(P ? P : 1);
72     switch (dim) {
73     case 1:
74       for (i = 0; i < isize+1; ++i) {
75         PetscInt v = i+vStart, dof, off;
76 
77         ierr = PetscSectionGetDof(csection, v, &dof);CHKERRQ(ierr);
78         ierr = PetscSectionGetOffset(csection, v, &off);CHKERRQ(ierr);
79         if (off >= 0) {
80           coors[off] = xmin + hx*(i+istart);
81         }
82       }
83       break;
84     case 2:
85       for (j = 0; j < jsize+1; ++j) {
86         for (i = 0; i < isize+1; ++i) {
87           PetscInt v = j*(isize+1)+i+vStart, dof, off;
88 
89           ierr = PetscSectionGetDof(csection, v, &dof);CHKERRQ(ierr);
90           ierr = PetscSectionGetOffset(csection, v, &off);CHKERRQ(ierr);
91           if (off >= 0) {
92             coors[off+0] = xmin + hx*(i+istart);
93             coors[off+1] = ymin + hy*(j+jstart);
94           }
95         }
96       }
97       break;
98     case 3:
99       for (k = 0; k < ksize+1; ++k) {
100         for (j = 0; j < jsize+1; ++j) {
101           for (i = 0; i < isize+1; ++i) {
102             PetscInt v = (k*(jsize+1)+j)*(isize+1)+i+vStart, dof, off;
103 
104             ierr = PetscSectionGetDof(csection, v, &dof);CHKERRQ(ierr);
105             ierr = PetscSectionGetOffset(csection, v, &off);CHKERRQ(ierr);
106             if (off >= 0) {
107               coors[off+0] = xmin + hx*(i+istart);
108               coors[off+1] = ymin + hy*(j+jstart);
109               coors[off+2] = zmin + hz_*(k+kstart);
110             }
111           }
112         }
113       }
114       break;
115     default:
116       SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Cannot create uniform coordinates for this dimension %D\n",dim);
117     }
118     ierr = VecRestoreArray(xcoor,&coors);CHKERRQ(ierr);
119     ierr = DMSetCoordinates(da,xcoor);CHKERRQ(ierr);
120     ierr = PetscLogObjectParent(da,xcoor);CHKERRQ(ierr);
121     ierr = VecDestroy(&xcoor);CHKERRQ(ierr);
122     PetscFunctionReturn(0);
123   }
124   if (dim == 1) {
125     if (bx == DMDA_BOUNDARY_PERIODIC) hx = (xmax-xmin)/M;
126     else                         hx = (xmax-xmin)/(M-1);
127     ierr = VecGetArray(xcoor,&coors);CHKERRQ(ierr);
128     for (i=0; i<isize; i++) {
129       coors[i] = xmin + hx*(i+istart);
130     }
131     ierr = VecRestoreArray(xcoor,&coors);CHKERRQ(ierr);
132   } else if (dim == 2) {
133     if (bx == DMDA_BOUNDARY_PERIODIC) hx = (xmax-xmin)/(M);
134     else                       hx = (xmax-xmin)/(M-1);
135     if (by == DMDA_BOUNDARY_PERIODIC) hy = (ymax-ymin)/(N);
136     else                       hy = (ymax-ymin)/(N-1);
137     ierr = VecGetArray(xcoor,&coors);CHKERRQ(ierr);
138     cnt  = 0;
139     for (j=0; j<jsize; j++) {
140       for (i=0; i<isize; i++) {
141         coors[cnt++] = xmin + hx*(i+istart);
142         coors[cnt++] = ymin + hy*(j+jstart);
143       }
144     }
145     ierr = VecRestoreArray(xcoor,&coors);CHKERRQ(ierr);
146   } else if (dim == 3) {
147     if (bx == DMDA_BOUNDARY_PERIODIC) hx = (xmax-xmin)/(M);
148     else                       hx = (xmax-xmin)/(M-1);
149     if (by == DMDA_BOUNDARY_PERIODIC) hy = (ymax-ymin)/(N);
150     else                       hy = (ymax-ymin)/(N-1);
151     if (bz == DMDA_BOUNDARY_PERIODIC) hz_ = (zmax-zmin)/(P);
152     else                       hz_ = (zmax-zmin)/(P-1);
153     ierr = VecGetArray(xcoor,&coors);CHKERRQ(ierr);
154     cnt  = 0;
155     for (k=0; k<ksize; k++) {
156       for (j=0; j<jsize; j++) {
157         for (i=0; i<isize; i++) {
158           coors[cnt++] = xmin + hx*(i+istart);
159           coors[cnt++] = ymin + hy*(j+jstart);
160           coors[cnt++] = zmin + hz_*(k+kstart);
161         }
162       }
163     }
164     ierr = VecRestoreArray(xcoor,&coors);CHKERRQ(ierr);
165   } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Cannot create uniform coordinates for this dimension %D\n",dim);
166   ierr = DMSetCoordinates(da,xcoor);CHKERRQ(ierr);
167   ierr = PetscLogObjectParent(da,xcoor);CHKERRQ(ierr);
168   ierr = VecDestroy(&xcoor);CHKERRQ(ierr);
169   PetscFunctionReturn(0);
170 }
171 
172 #undef __FUNCT__
173 #define __FUNCT__ "VecView_MPI_Draw_DA1d"
174 PetscErrorCode VecView_MPI_Draw_DA1d(Vec xin,PetscViewer v)
175 {
176   DM                da;
177   PetscErrorCode    ierr;
178   PetscMPIInt       rank,size,tag1,tag2;
179   PetscInt          i,n,N,step,istart,isize,j,nbounds;
180   MPI_Status        status;
181   PetscReal         coors[4],ymin,ymax,min,max,xmin,xmax,tmp,xgtmp;
182   const PetscScalar *array,*xg;
183   PetscDraw         draw;
184   PetscBool         isnull,showpoints = PETSC_FALSE;
185   MPI_Comm          comm;
186   PetscDrawAxis     axis;
187   Vec               xcoor;
188   DMDABoundaryType  bx;
189   const PetscReal   *bounds;
190   PetscInt          *displayfields;
191   PetscInt          k,ndisplayfields;
192   PetscBool         flg,hold;
193 
194   PetscFunctionBegin;
195   ierr = PetscViewerDrawGetDraw(v,0,&draw);CHKERRQ(ierr);
196   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
197   ierr = PetscViewerDrawGetBounds(v,&nbounds,&bounds);CHKERRQ(ierr);
198 
199   ierr = PetscObjectQuery((PetscObject)xin,"DM",(PetscObject*)&da);CHKERRQ(ierr);
200   if (!da) SETERRQ(((PetscObject)xin)->comm,PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
201 
202   ierr = PetscOptionsGetBool(PETSC_NULL,"-draw_vec_mark_points",&showpoints,PETSC_NULL);CHKERRQ(ierr);
203 
204   ierr = DMDAGetInfo(da,0,&N,0,0,0,0,0,&step,0,&bx,0,0,0);CHKERRQ(ierr);
205   ierr = DMDAGetCorners(da,&istart,0,0,&isize,0,0);CHKERRQ(ierr);
206   ierr = VecGetArrayRead(xin,&array);CHKERRQ(ierr);
207   ierr = VecGetLocalSize(xin,&n);CHKERRQ(ierr);
208   n    = n/step;
209 
210   /* get coordinates of nodes */
211   ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr);
212   if (!xcoor) {
213     ierr = DMDASetUniformCoordinates(da,0.0,1.0,0.0,0.0,0.0,0.0);CHKERRQ(ierr);
214     ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr);
215   }
216   ierr = VecGetArrayRead(xcoor,&xg);CHKERRQ(ierr);
217 
218   ierr = PetscObjectGetComm((PetscObject)xin,&comm);CHKERRQ(ierr);
219   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
220   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
221 
222   /*
223       Determine the min and max x coordinate in plot
224   */
225   if (!rank) {
226     xmin = PetscRealPart(xg[0]);
227   }
228   if (rank == size-1) {
229     xmax = PetscRealPart(xg[n-1]);
230   }
231   ierr = MPI_Bcast(&xmin,1,MPIU_REAL,0,comm);CHKERRQ(ierr);
232   ierr = MPI_Bcast(&xmax,1,MPIU_REAL,size-1,comm);CHKERRQ(ierr);
233 
234   ierr = PetscMalloc(step*sizeof(PetscInt),&displayfields);CHKERRQ(ierr);
235   for (i=0; i<step; i++) displayfields[i] = i;
236   ndisplayfields = step;
237   ierr = PetscOptionsGetIntArray(PETSC_NULL,"-draw_fields",displayfields,&ndisplayfields,&flg);CHKERRQ(ierr);
238   if (!flg) ndisplayfields = step;
239   for (k=0; k<ndisplayfields; k++) {
240     j = displayfields[k];
241     ierr = PetscViewerDrawGetDraw(v,k,&draw);CHKERRQ(ierr);
242     ierr = PetscDrawCheckResizedWindow(draw);CHKERRQ(ierr);
243 
244     /*
245         Determine the min and max y coordinate in plot
246     */
247     min = 1.e20; max = -1.e20;
248     for (i=0; i<n; i++) {
249       if (PetscRealPart(array[j+i*step]) < min) min = PetscRealPart(array[j+i*step]);
250       if (PetscRealPart(array[j+i*step]) > max) max = PetscRealPart(array[j+i*step]);
251     }
252     if (min + 1.e-10 > max) {
253       min -= 1.e-5;
254       max += 1.e-5;
255     }
256     if (j < nbounds) {
257       min = PetscMin(min,bounds[2*j]);
258       max = PetscMax(max,bounds[2*j+1]);
259     }
260 
261     ierr = MPI_Reduce(&min,&ymin,1,MPIU_REAL,MPIU_MIN,0,comm);CHKERRQ(ierr);
262     ierr = MPI_Reduce(&max,&ymax,1,MPIU_REAL,MPIU_MAX,0,comm);CHKERRQ(ierr);
263 
264     ierr = PetscViewerDrawGetHold(v,&hold);CHKERRQ(ierr);
265     if (!hold) {
266       ierr = PetscDrawSynchronizedClear(draw);CHKERRQ(ierr);
267     }
268     ierr = PetscViewerDrawGetDrawAxis(v,k,&axis);CHKERRQ(ierr);
269     ierr = PetscLogObjectParent(draw,axis);CHKERRQ(ierr);
270     if (!rank) {
271       const char *title;
272 
273       ierr = PetscDrawAxisSetLimits(axis,xmin,xmax,ymin,ymax);CHKERRQ(ierr);
274       ierr = PetscDrawAxisDraw(axis);CHKERRQ(ierr);
275       ierr = PetscDrawGetCoordinates(draw,coors,coors+1,coors+2,coors+3);CHKERRQ(ierr);
276       ierr = DMDAGetFieldName(da,j,&title);CHKERRQ(ierr);
277       if (title) {ierr = PetscDrawSetTitle(draw,title);CHKERRQ(ierr);}
278     }
279     ierr = MPI_Bcast(coors,4,MPIU_REAL,0,comm);CHKERRQ(ierr);
280     if (rank) {
281       ierr = PetscDrawSetCoordinates(draw,coors[0],coors[1],coors[2],coors[3]);CHKERRQ(ierr);
282     }
283 
284     /* draw local part of vector */
285     ierr = PetscObjectGetNewTag((PetscObject)xin,&tag1);CHKERRQ(ierr);
286     ierr = PetscObjectGetNewTag((PetscObject)xin,&tag2);CHKERRQ(ierr);
287     if (rank < size-1) { /*send value to right */
288       ierr = MPI_Send((void*)&array[j+(n-1)*step],1,MPIU_REAL,rank+1,tag1,comm);CHKERRQ(ierr);
289       ierr = MPI_Send((void*)&xg[n-1],1,MPIU_REAL,rank+1,tag1,comm);CHKERRQ(ierr);
290     }
291     if (!rank && bx == DMDA_BOUNDARY_PERIODIC && size > 1) { /* first processor sends first value to last */
292       ierr = MPI_Send((void*)&array[j],1,MPIU_REAL,size-1,tag2,comm);CHKERRQ(ierr);
293     }
294 
295     for (i=1; i<n; i++) {
296       ierr = PetscDrawLine(draw,PetscRealPart(xg[i-1]),PetscRealPart(array[j+step*(i-1)]),PetscRealPart(xg[i]),PetscRealPart(array[j+step*i]),PETSC_DRAW_RED);CHKERRQ(ierr);
297       if (showpoints) {
298         ierr = PetscDrawPoint(draw,PetscRealPart(xg[i-1]),PetscRealPart(array[j+step*(i-1)]),PETSC_DRAW_BLACK);CHKERRQ(ierr);
299       }
300     }
301     if (rank) { /* receive value from left */
302       ierr = MPI_Recv(&tmp,1,MPIU_REAL,rank-1,tag1,comm,&status);CHKERRQ(ierr);
303       ierr = MPI_Recv(&xgtmp,1,MPIU_REAL,rank-1,tag1,comm,&status);CHKERRQ(ierr);
304       ierr = PetscDrawLine(draw,xgtmp,tmp,PetscRealPart(xg[0]),PetscRealPart(array[j]),PETSC_DRAW_RED);CHKERRQ(ierr);
305       if (showpoints) {
306         ierr = PetscDrawPoint(draw,xgtmp,tmp,PETSC_DRAW_BLACK);CHKERRQ(ierr);
307       }
308     }
309     if (rank == size-1 && bx == DMDA_BOUNDARY_PERIODIC && size > 1) {
310       ierr = MPI_Recv(&tmp,1,MPIU_REAL,0,tag2,comm,&status);CHKERRQ(ierr);
311       /* If the mesh is not uniform we do not know the mesh spacing between the last point on the right and the first ghost point */
312       ierr = PetscDrawLine(draw,PetscRealPart(xg[n-1]),PetscRealPart(array[j+step*(n-1)]),PetscRealPart(xg[n-1]+(xg[n-1]-xg[n-2])),tmp,PETSC_DRAW_RED);CHKERRQ(ierr);
313       if (showpoints) {
314         ierr = PetscDrawPoint(draw,PetscRealPart(xg[n-2]),PetscRealPart(array[j+step*(n-1)]),PETSC_DRAW_BLACK);CHKERRQ(ierr);
315       }
316     }
317     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
318     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
319   }
320   ierr = PetscFree(displayfields);CHKERRQ(ierr);
321   ierr = VecRestoreArrayRead(xcoor,&xg);CHKERRQ(ierr);
322   ierr = VecRestoreArrayRead(xin,&array);CHKERRQ(ierr);
323   PetscFunctionReturn(0);
324 }
325 
326