xref: /petsc/src/dm/impls/da/gr1.c (revision 750b007cd8d816cecd9de99077bb0a703b4cf61a)
1 
2 /*
3    Plots vectors obtained with DMDACreate1d()
4 */
5 
6 #include <petsc/private/dmdaimpl.h> /*I  "petscdmda.h"   I*/
7 
8 /*@
9     DMDASetUniformCoordinates - Sets a DMDA coordinates to be a uniform grid
10 
11   Collective on da
12 
13   Input Parameters:
14 +  da - the distributed array object
15 .  xmin,xmax - extremes in the x direction
16 .  ymin,ymax - extremes in the y direction (value ignored for 1 dimensional problems)
17 -  zmin,zmax - extremes in the z direction (value ignored for 1 or 2 dimensional problems)
18 
19   Level: beginner
20 
21 .seealso: `DMSetCoordinates()`, `DMGetCoordinates()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`, `DMStagSetUniformCoordinates()`
22 
23 @*/
24 PetscErrorCode DMDASetUniformCoordinates(DM da, PetscReal xmin, PetscReal xmax, PetscReal ymin, PetscReal ymax, PetscReal zmin, PetscReal zmax) {
25   MPI_Comm       comm;
26   DM             cda;
27   DM_DA         *dd = (DM_DA *)da->data;
28   DMBoundaryType bx, by, bz;
29   Vec            xcoor;
30   PetscScalar   *coors;
31   PetscReal      hx, hy, hz_;
32   PetscInt       i, j, k, M, N, P, istart, isize, jstart, jsize, kstart, ksize, dim, cnt;
33 
34   PetscFunctionBegin;
35   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
36   PetscCheck(dd->gtol, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "Cannot set coordinates until after DMDA has been setup");
37   PetscCall(DMDAGetInfo(da, &dim, &M, &N, &P, NULL, NULL, NULL, NULL, NULL, &bx, &by, &bz, NULL));
38   PetscCheck(xmax >= xmin, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_INCOMP, "xmax must be larger than xmin %g %g", (double)xmin, (double)xmax);
39   PetscCheck(!(dim > 1) || !(ymax < ymin), PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_INCOMP, "ymax must be larger than ymin %g %g", (double)ymin, (double)ymax);
40   PetscCheck(!(dim > 2) || !(zmax < zmin), PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_INCOMP, "zmax must be larger than zmin %g %g", (double)zmin, (double)zmax);
41   PetscCall(PetscObjectGetComm((PetscObject)da, &comm));
42   PetscCall(DMDAGetCorners(da, &istart, &jstart, &kstart, &isize, &jsize, &ksize));
43   PetscCall(DMGetCoordinateDM(da, &cda));
44   PetscCall(DMCreateGlobalVector(cda, &xcoor));
45   if (dim == 1) {
46     if (bx == DM_BOUNDARY_PERIODIC) hx = (xmax - xmin) / M;
47     else hx = (xmax - xmin) / (M - 1);
48     PetscCall(VecGetArray(xcoor, &coors));
49     for (i = 0; i < isize; i++) coors[i] = xmin + hx * (i + istart);
50     PetscCall(VecRestoreArray(xcoor, &coors));
51   } else if (dim == 2) {
52     if (bx == DM_BOUNDARY_PERIODIC) hx = (xmax - xmin) / (M);
53     else hx = (xmax - xmin) / (M - 1);
54     if (by == DM_BOUNDARY_PERIODIC) hy = (ymax - ymin) / (N);
55     else hy = (ymax - ymin) / (N - 1);
56     PetscCall(VecGetArray(xcoor, &coors));
57     cnt = 0;
58     for (j = 0; j < jsize; j++) {
59       for (i = 0; i < isize; i++) {
60         coors[cnt++] = xmin + hx * (i + istart);
61         coors[cnt++] = ymin + hy * (j + jstart);
62       }
63     }
64     PetscCall(VecRestoreArray(xcoor, &coors));
65   } else if (dim == 3) {
66     if (bx == DM_BOUNDARY_PERIODIC) hx = (xmax - xmin) / (M);
67     else hx = (xmax - xmin) / (M - 1);
68     if (by == DM_BOUNDARY_PERIODIC) hy = (ymax - ymin) / (N);
69     else hy = (ymax - ymin) / (N - 1);
70     if (bz == DM_BOUNDARY_PERIODIC) hz_ = (zmax - zmin) / (P);
71     else hz_ = (zmax - zmin) / (P - 1);
72     PetscCall(VecGetArray(xcoor, &coors));
73     cnt = 0;
74     for (k = 0; k < ksize; k++) {
75       for (j = 0; j < jsize; j++) {
76         for (i = 0; i < isize; i++) {
77           coors[cnt++] = xmin + hx * (i + istart);
78           coors[cnt++] = ymin + hy * (j + jstart);
79           coors[cnt++] = zmin + hz_ * (k + kstart);
80         }
81       }
82     }
83     PetscCall(VecRestoreArray(xcoor, &coors));
84   } else SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Cannot create uniform coordinates for this dimension %" PetscInt_FMT, dim);
85   PetscCall(DMSetCoordinates(da, xcoor));
86   PetscCall(VecDestroy(&xcoor));
87   PetscFunctionReturn(0);
88 }
89 
90 /*
91     Allows a user to select a subset of the fields to be drawn by VecView() when the vector comes from a DMDA
92 */
93 PetscErrorCode DMDASelectFields(DM da, PetscInt *outfields, PetscInt **fields) {
94   PetscInt  step, ndisplayfields, *displayfields, k, j;
95   PetscBool flg;
96 
97   PetscFunctionBegin;
98   PetscCall(DMDAGetInfo(da, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &step, NULL, NULL, NULL, NULL, NULL));
99   PetscCall(PetscMalloc1(step, &displayfields));
100   for (k = 0; k < step; k++) displayfields[k] = k;
101   ndisplayfields = step;
102   PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-draw_fields", displayfields, &ndisplayfields, &flg));
103   if (!ndisplayfields) ndisplayfields = step;
104   if (!flg) {
105     char      **fields;
106     const char *fieldname;
107     PetscInt    nfields = step;
108     PetscCall(PetscMalloc1(step, &fields));
109     PetscCall(PetscOptionsGetStringArray(NULL, NULL, "-draw_fields_by_name", fields, &nfields, &flg));
110     if (flg) {
111       ndisplayfields = 0;
112       for (k = 0; k < nfields; k++) {
113         for (j = 0; j < step; j++) {
114           PetscCall(DMDAGetFieldName(da, j, &fieldname));
115           PetscCall(PetscStrcmp(fieldname, fields[k], &flg));
116           if (flg) goto found;
117         }
118         SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_USER, "Unknown fieldname %s", fields[k]);
119       found:
120         displayfields[ndisplayfields++] = j;
121       }
122     }
123     for (k = 0; k < nfields; k++) PetscCall(PetscFree(fields[k]));
124     PetscCall(PetscFree(fields));
125   }
126   *fields    = displayfields;
127   *outfields = ndisplayfields;
128   PetscFunctionReturn(0);
129 }
130 
131 #include <petscdraw.h>
132 
133 PetscErrorCode VecView_MPI_Draw_DA1d(Vec xin, PetscViewer v) {
134   DM                  da;
135   PetscMPIInt         rank, size, tag;
136   PetscInt            i, n, N, dof, istart, isize, j, nbounds;
137   MPI_Status          status;
138   PetscReal           min, max, xmin = 0.0, xmax = 0.0, tmp = 0.0, xgtmp = 0.0;
139   const PetscScalar  *array, *xg;
140   PetscDraw           draw;
141   PetscBool           isnull, useports = PETSC_FALSE, showmarkers = PETSC_FALSE;
142   MPI_Comm            comm;
143   PetscDrawAxis       axis;
144   Vec                 xcoor;
145   DMBoundaryType      bx;
146   const char         *tlabel = NULL, *xlabel = NULL;
147   const PetscReal    *bounds;
148   PetscInt           *displayfields;
149   PetscInt            k, ndisplayfields;
150   PetscBool           hold;
151   PetscDrawViewPorts *ports = NULL;
152   PetscViewerFormat   format;
153 
154   PetscFunctionBegin;
155   PetscCall(PetscViewerDrawGetDraw(v, 0, &draw));
156   PetscCall(PetscDrawIsNull(draw, &isnull));
157   if (isnull) PetscFunctionReturn(0);
158   PetscCall(PetscViewerDrawGetBounds(v, &nbounds, &bounds));
159 
160   PetscCall(VecGetDM(xin, &da));
161   PetscCheck(da, PetscObjectComm((PetscObject)xin), PETSC_ERR_ARG_WRONG, "Vector not generated from a DMDA");
162   PetscCall(PetscObjectGetComm((PetscObject)xin, &comm));
163   PetscCallMPI(MPI_Comm_size(comm, &size));
164   PetscCallMPI(MPI_Comm_rank(comm, &rank));
165 
166   PetscCall(PetscOptionsGetBool(NULL, NULL, "-draw_vec_use_markers", &showmarkers, NULL));
167 
168   PetscCall(DMDAGetInfo(da, NULL, &N, NULL, NULL, NULL, NULL, NULL, &dof, NULL, &bx, NULL, NULL, NULL));
169   PetscCall(DMDAGetCorners(da, &istart, NULL, NULL, &isize, NULL, NULL));
170   PetscCall(VecGetArrayRead(xin, &array));
171   PetscCall(VecGetLocalSize(xin, &n));
172   n = n / dof;
173 
174   /* Get coordinates of nodes */
175   PetscCall(DMGetCoordinates(da, &xcoor));
176   if (!xcoor) {
177     PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 0.0, 0.0, 0.0));
178     PetscCall(DMGetCoordinates(da, &xcoor));
179   }
180   PetscCall(VecGetArrayRead(xcoor, &xg));
181   PetscCall(DMDAGetCoordinateName(da, 0, &xlabel));
182 
183   /* Determine the min and max coordinate in plot */
184   if (rank == 0) xmin = PetscRealPart(xg[0]);
185   if (rank == size - 1) xmax = PetscRealPart(xg[n - 1]);
186   PetscCallMPI(MPI_Bcast(&xmin, 1, MPIU_REAL, 0, comm));
187   PetscCallMPI(MPI_Bcast(&xmax, 1, MPIU_REAL, size - 1, comm));
188 
189   PetscCall(DMDASelectFields(da, &ndisplayfields, &displayfields));
190   PetscCall(PetscViewerGetFormat(v, &format));
191   PetscCall(PetscOptionsGetBool(NULL, NULL, "-draw_ports", &useports, NULL));
192   if (format == PETSC_VIEWER_DRAW_PORTS) useports = PETSC_TRUE;
193   if (useports) {
194     PetscCall(PetscViewerDrawGetDraw(v, 0, &draw));
195     PetscCall(PetscViewerDrawGetDrawAxis(v, 0, &axis));
196     PetscCall(PetscDrawCheckResizedWindow(draw));
197     PetscCall(PetscDrawClear(draw));
198     PetscCall(PetscDrawViewPortsCreate(draw, ndisplayfields, &ports));
199   }
200 
201   /* Loop over each field; drawing each in a different window */
202   for (k = 0; k < ndisplayfields; k++) {
203     j = displayfields[k];
204 
205     /* determine the min and max value in plot */
206     PetscCall(VecStrideMin(xin, j, NULL, &min));
207     PetscCall(VecStrideMax(xin, j, NULL, &max));
208     if (j < nbounds) {
209       min = PetscMin(min, bounds[2 * j]);
210       max = PetscMax(max, bounds[2 * j + 1]);
211     }
212     if (min == max) {
213       min -= 1.e-5;
214       max += 1.e-5;
215     }
216 
217     if (useports) {
218       PetscCall(PetscDrawViewPortsSet(ports, k));
219       PetscCall(DMDAGetFieldName(da, j, &tlabel));
220     } else {
221       const char *title;
222       PetscCall(PetscViewerDrawGetHold(v, &hold));
223       PetscCall(PetscViewerDrawGetDraw(v, k, &draw));
224       PetscCall(PetscViewerDrawGetDrawAxis(v, k, &axis));
225       PetscCall(DMDAGetFieldName(da, j, &title));
226       if (title) PetscCall(PetscDrawSetTitle(draw, title));
227       PetscCall(PetscDrawCheckResizedWindow(draw));
228       if (!hold) PetscCall(PetscDrawClear(draw));
229     }
230     PetscCall(PetscDrawAxisSetLabels(axis, tlabel, xlabel, NULL));
231     PetscCall(PetscDrawAxisSetLimits(axis, xmin, xmax, min, max));
232     PetscCall(PetscDrawAxisDraw(axis));
233 
234     /* draw local part of vector */
235     PetscCall(PetscObjectGetNewTag((PetscObject)xin, &tag));
236     if (rank < size - 1) { /*send value to right */
237       PetscCallMPI(MPI_Send((void *)&xg[n - 1], 1, MPIU_REAL, rank + 1, tag, comm));
238       PetscCallMPI(MPI_Send((void *)&array[j + (n - 1) * dof], 1, MPIU_REAL, rank + 1, tag, comm));
239     }
240     if (rank) { /* receive value from left */
241       PetscCallMPI(MPI_Recv(&xgtmp, 1, MPIU_REAL, rank - 1, tag, comm, &status));
242       PetscCallMPI(MPI_Recv(&tmp, 1, MPIU_REAL, rank - 1, tag, comm, &status));
243     }
244     PetscDrawCollectiveBegin(draw);
245     if (rank) {
246       PetscCall(PetscDrawLine(draw, xgtmp, tmp, PetscRealPart(xg[0]), PetscRealPart(array[j]), PETSC_DRAW_RED));
247       if (showmarkers) PetscCall(PetscDrawPoint(draw, xgtmp, tmp, PETSC_DRAW_BLACK));
248     }
249     for (i = 1; i < n; i++) {
250       PetscCall(PetscDrawLine(draw, PetscRealPart(xg[i - 1]), PetscRealPart(array[j + dof * (i - 1)]), PetscRealPart(xg[i]), PetscRealPart(array[j + dof * i]), PETSC_DRAW_RED));
251       if (showmarkers) PetscCall(PetscDrawMarker(draw, PetscRealPart(xg[i - 1]), PetscRealPart(array[j + dof * (i - 1)]), PETSC_DRAW_BLACK));
252     }
253     if (rank == size - 1) {
254       if (showmarkers) PetscCall(PetscDrawMarker(draw, PetscRealPart(xg[n - 1]), PetscRealPart(array[j + dof * (n - 1)]), PETSC_DRAW_BLACK));
255     }
256     PetscDrawCollectiveEnd(draw);
257     PetscCall(PetscDrawFlush(draw));
258     PetscCall(PetscDrawPause(draw));
259     if (!useports) PetscCall(PetscDrawSave(draw));
260   }
261   if (useports) {
262     PetscCall(PetscViewerDrawGetDraw(v, 0, &draw));
263     PetscCall(PetscDrawSave(draw));
264   }
265 
266   PetscCall(PetscDrawViewPortsDestroy(ports));
267   PetscCall(PetscFree(displayfields));
268   PetscCall(VecRestoreArrayRead(xcoor, &xg));
269   PetscCall(VecRestoreArrayRead(xin, &array));
270   PetscFunctionReturn(0);
271 }
272