xref: /petsc/src/dm/impls/da/gr1.c (revision 2205254efee3a00a594e5e2a3a70f74dcb40bc03)
1 
2 /*
3    Plots vectors obtained with DMDACreate1d()
4 */
5 
6 #include <petsc-private/daimpl.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__ "DMDASelectFields"
174 PetscErrorCode DMDASelectFields(DM da,PetscInt *outfields,PetscInt **fields)
175 {
176   PetscErrorCode ierr;
177   PetscInt       step,ndisplayfields,*displayfields,k,j;
178   PetscBool      flg;
179 
180   PetscFunctionBegin;
181   ierr = DMDAGetInfo(da,0,0,0,0,0,0,0,&step,0,0,0,0,0);CHKERRQ(ierr);
182   ierr = PetscMalloc(step*sizeof(PetscInt),&displayfields);CHKERRQ(ierr);
183   for (k=0; k<step; k++) displayfields[k] = k;
184   ndisplayfields = step;
185   ierr           = PetscOptionsGetIntArray(PETSC_NULL,"-draw_fields",displayfields,&ndisplayfields,&flg);CHKERRQ(ierr);
186   if (!ndisplayfields) ndisplayfields = step;
187   if (!flg) {
188     char       **fields;
189     const char *fieldname;
190     PetscInt   nfields = step;
191     ierr = PetscMalloc(step*sizeof(char*),&fields);CHKERRQ(ierr);
192     ierr = PetscOptionsGetStringArray(PETSC_NULL,"-draw_fields_by_name",fields,&nfields,&flg);CHKERRQ(ierr);
193     if (flg) {
194       ndisplayfields = 0;
195       for (k=0; k<nfields;k++) {
196         for (j=0; j<step; j++) {
197           ierr = DMDAGetFieldName(da,j,&fieldname);CHKERRQ(ierr);
198           ierr = PetscStrcmp(fieldname,fields[k],&flg);CHKERRQ(ierr);
199           if (flg) {
200             goto found;
201           }
202         }
203         SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_USER,"Unknown fieldname %s",fields[k]);
204 found:  displayfields[ndisplayfields++] = j;
205       }
206     }
207     for (k=0; k<nfields; k++) {
208       ierr = PetscFree(fields[k]);CHKERRQ(ierr);
209     }
210     ierr = PetscFree(fields);CHKERRQ(ierr);
211   }
212   *fields    = displayfields;
213   *outfields = ndisplayfields;
214   PetscFunctionReturn(0);
215 }
216 
217 #undef __FUNCT__
218 #define __FUNCT__ "VecView_MPI_Draw_DA1d"
219 PetscErrorCode VecView_MPI_Draw_DA1d(Vec xin,PetscViewer v)
220 {
221   DM                da;
222   PetscErrorCode    ierr;
223   PetscMPIInt       rank,size,tag1,tag2;
224   PetscInt          i,n,N,step,istart,isize,j,nbounds;
225   MPI_Status        status;
226   PetscReal         coors[4],ymin,ymax,min,max,xmin,xmax,tmp,xgtmp;
227   const PetscScalar *array,*xg;
228   PetscDraw         draw;
229   PetscBool         isnull,showpoints = PETSC_FALSE;
230   MPI_Comm          comm;
231   PetscDrawAxis     axis;
232   Vec               xcoor;
233   DMDABoundaryType  bx;
234   const PetscReal   *bounds;
235   PetscInt          *displayfields;
236   PetscInt          k,ndisplayfields;
237   PetscBool         hold;
238 
239   PetscFunctionBegin;
240   ierr = PetscViewerDrawGetDraw(v,0,&draw);CHKERRQ(ierr);
241   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
242   ierr = PetscViewerDrawGetBounds(v,&nbounds,&bounds);CHKERRQ(ierr);
243 
244   ierr = VecGetDM(xin,&da);CHKERRQ(ierr);
245   if (!da) SETERRQ(((PetscObject)xin)->comm,PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
246 
247   ierr = PetscOptionsGetBool(PETSC_NULL,"-draw_vec_mark_points",&showpoints,PETSC_NULL);CHKERRQ(ierr);
248 
249   ierr = DMDAGetInfo(da,0,&N,0,0,0,0,0,&step,0,&bx,0,0,0);CHKERRQ(ierr);
250   ierr = DMDAGetCorners(da,&istart,0,0,&isize,0,0);CHKERRQ(ierr);
251   ierr = VecGetArrayRead(xin,&array);CHKERRQ(ierr);
252   ierr = VecGetLocalSize(xin,&n);CHKERRQ(ierr);
253   n    = n/step;
254 
255   /* get coordinates of nodes */
256   ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr);
257   if (!xcoor) {
258     ierr = DMDASetUniformCoordinates(da,0.0,1.0,0.0,0.0,0.0,0.0);CHKERRQ(ierr);
259     ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr);
260   }
261   ierr = VecGetArrayRead(xcoor,&xg);CHKERRQ(ierr);
262 
263   ierr = PetscObjectGetComm((PetscObject)xin,&comm);CHKERRQ(ierr);
264   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
265   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
266 
267   /*
268       Determine the min and max x coordinate in plot
269   */
270   if (!rank) {
271     xmin = PetscRealPart(xg[0]);
272   }
273   if (rank == size-1) {
274     xmax = PetscRealPart(xg[n-1]);
275   }
276   ierr = MPI_Bcast(&xmin,1,MPIU_REAL,0,comm);CHKERRQ(ierr);
277   ierr = MPI_Bcast(&xmax,1,MPIU_REAL,size-1,comm);CHKERRQ(ierr);
278 
279   ierr = DMDASelectFields(da,&ndisplayfields,&displayfields);CHKERRQ(ierr);
280   for (k=0; k<ndisplayfields; k++) {
281     j    = displayfields[k];
282     ierr = PetscViewerDrawGetDraw(v,k,&draw);CHKERRQ(ierr);
283     ierr = PetscDrawCheckResizedWindow(draw);CHKERRQ(ierr);
284 
285     /*
286         Determine the min and max y coordinate in plot
287     */
288     min = 1.e20; max = -1.e20;
289     for (i=0; i<n; i++) {
290       if (PetscRealPart(array[j+i*step]) < min) min = PetscRealPart(array[j+i*step]);
291       if (PetscRealPart(array[j+i*step]) > max) max = PetscRealPart(array[j+i*step]);
292     }
293     if (min + 1.e-10 > max) {
294       min -= 1.e-5;
295       max += 1.e-5;
296     }
297     if (j < nbounds) {
298       min = PetscMin(min,bounds[2*j]);
299       max = PetscMax(max,bounds[2*j+1]);
300     }
301 
302     ierr = MPI_Reduce(&min,&ymin,1,MPIU_REAL,MPIU_MIN,0,comm);CHKERRQ(ierr);
303     ierr = MPI_Reduce(&max,&ymax,1,MPIU_REAL,MPIU_MAX,0,comm);CHKERRQ(ierr);
304 
305     ierr = PetscViewerDrawGetHold(v,&hold);CHKERRQ(ierr);
306     if (!hold) {
307       ierr = PetscDrawSynchronizedClear(draw);CHKERRQ(ierr);
308     }
309     ierr = PetscViewerDrawGetDrawAxis(v,k,&axis);CHKERRQ(ierr);
310     ierr = PetscLogObjectParent(draw,axis);CHKERRQ(ierr);
311     if (!rank) {
312       const char *title;
313 
314       ierr = PetscDrawAxisSetLimits(axis,xmin,xmax,ymin,ymax);CHKERRQ(ierr);
315       ierr = PetscDrawAxisDraw(axis);CHKERRQ(ierr);
316       ierr = PetscDrawGetCoordinates(draw,coors,coors+1,coors+2,coors+3);CHKERRQ(ierr);
317       ierr = DMDAGetFieldName(da,j,&title);CHKERRQ(ierr);
318       if (title) {ierr = PetscDrawSetTitle(draw,title);CHKERRQ(ierr);}
319     }
320     ierr = MPI_Bcast(coors,4,MPIU_REAL,0,comm);CHKERRQ(ierr);
321     if (rank) {
322       ierr = PetscDrawSetCoordinates(draw,coors[0],coors[1],coors[2],coors[3]);CHKERRQ(ierr);
323     }
324 
325     /* draw local part of vector */
326     ierr = PetscObjectGetNewTag((PetscObject)xin,&tag1);CHKERRQ(ierr);
327     ierr = PetscObjectGetNewTag((PetscObject)xin,&tag2);CHKERRQ(ierr);
328     if (rank < size-1) { /*send value to right */
329       ierr = MPI_Send((void*)&array[j+(n-1)*step],1,MPIU_REAL,rank+1,tag1,comm);CHKERRQ(ierr);
330       ierr = MPI_Send((void*)&xg[n-1],1,MPIU_REAL,rank+1,tag1,comm);CHKERRQ(ierr);
331     }
332     if (!rank && bx == DMDA_BOUNDARY_PERIODIC && size > 1) { /* first processor sends first value to last */
333       ierr = MPI_Send((void*)&array[j],1,MPIU_REAL,size-1,tag2,comm);CHKERRQ(ierr);
334     }
335 
336     for (i=1; i<n; i++) {
337       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);
338       if (showpoints) {
339         ierr = PetscDrawPoint(draw,PetscRealPart(xg[i-1]),PetscRealPart(array[j+step*(i-1)]),PETSC_DRAW_BLACK);CHKERRQ(ierr);
340       }
341     }
342     if (rank) { /* receive value from left */
343       ierr = MPI_Recv(&tmp,1,MPIU_REAL,rank-1,tag1,comm,&status);CHKERRQ(ierr);
344       ierr = MPI_Recv(&xgtmp,1,MPIU_REAL,rank-1,tag1,comm,&status);CHKERRQ(ierr);
345       ierr = PetscDrawLine(draw,xgtmp,tmp,PetscRealPart(xg[0]),PetscRealPart(array[j]),PETSC_DRAW_RED);CHKERRQ(ierr);
346       if (showpoints) {
347         ierr = PetscDrawPoint(draw,xgtmp,tmp,PETSC_DRAW_BLACK);CHKERRQ(ierr);
348       }
349     }
350     if (rank == size-1 && bx == DMDA_BOUNDARY_PERIODIC && size > 1) {
351       ierr = MPI_Recv(&tmp,1,MPIU_REAL,0,tag2,comm,&status);CHKERRQ(ierr);
352       /* 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 */
353       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);
354       if (showpoints) {
355         ierr = PetscDrawPoint(draw,PetscRealPart(xg[n-2]),PetscRealPart(array[j+step*(n-1)]),PETSC_DRAW_BLACK);CHKERRQ(ierr);
356       }
357     }
358     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
359     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
360   }
361   ierr = PetscFree(displayfields);CHKERRQ(ierr);
362   ierr = VecRestoreArrayRead(xcoor,&xg);CHKERRQ(ierr);
363   ierr = VecRestoreArrayRead(xin,&array);CHKERRQ(ierr);
364   PetscFunctionReturn(0);
365 }
366 
367