Actual source code: gr1.c

  1: /*
  2:    Plots vectors obtained with DMDACreate1d()
  3: */

  5: #include <petsc/private/dmdaimpl.h>

  7: /*@
  8:   DMDASetUniformCoordinates - Sets a `DMDA` coordinates to be a uniform grid

 10:   Collective

 12:   Input Parameters:
 13: + da   - the distributed array object
 14: . xmin - min extreme in the x direction
 15: . xmax - max extreme in the x direction
 16: . ymin - min extreme in the y direction (value ignored for 1 dimensional problems)
 17: . ymax - max extreme in the y direction (value ignored for 1 dimensional problems)
 18: . zmin - min extreme in the z direction (value ignored for 1 or 2 dimensional problems)
 19: - zmax - max extreme in the z direction (value ignored for 1 or 2 dimensional problems)

 21:   Level: beginner

 23: .seealso: [](sec_struct), `DM`, `DMDA`, `DMSetCoordinates()`, `DMGetCoordinates()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`, `DMStagSetUniformCoordinates()`
 24: @*/
 25: PetscErrorCode DMDASetUniformCoordinates(DM da, PetscReal xmin, PetscReal xmax, PetscReal ymin, PetscReal ymax, PetscReal zmin, PetscReal zmax)
 26: {
 27:   MPI_Comm       comm;
 28:   DM             cda;
 29:   DM_DA         *dd = (DM_DA *)da->data;
 30:   DMBoundaryType bx, by, bz;
 31:   Vec            xcoor;
 32:   PetscScalar   *coors;
 33:   PetscReal      hx, hy, hz_;
 34:   PetscInt       i, j, k, M, N, P, istart, isize, jstart, jsize, kstart, ksize, dim, cnt;

 36:   PetscFunctionBegin;
 38:   PetscCheck(dd->gtol, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "Cannot set coordinates until after DMDA has been setup");
 39:   PetscCall(DMDAGetInfo(da, &dim, &M, &N, &P, NULL, NULL, NULL, NULL, NULL, &bx, &by, &bz, NULL));
 40:   PetscCheck(xmax >= xmin, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_INCOMP, "xmax must be larger than xmin %g %g", (double)xmin, (double)xmax);
 41:   PetscCheck(!(dim > 1) || !(ymax < ymin), PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_INCOMP, "ymax must be larger than ymin %g %g", (double)ymin, (double)ymax);
 42:   PetscCheck(!(dim > 2) || !(zmax < zmin), PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_INCOMP, "zmax must be larger than zmin %g %g", (double)zmin, (double)zmax);
 43:   PetscCall(PetscObjectGetComm((PetscObject)da, &comm));
 44:   PetscCall(DMDAGetCorners(da, &istart, &jstart, &kstart, &isize, &jsize, &ksize));
 45:   PetscCall(DMGetCoordinateDM(da, &cda));
 46:   PetscCall(DMCreateGlobalVector(cda, &xcoor));
 47:   if (dim == 1) {
 48:     if (bx == DM_BOUNDARY_PERIODIC) hx = (xmax - xmin) / M;
 49:     else hx = (xmax - xmin) / (M - 1);
 50:     PetscCall(VecGetArray(xcoor, &coors));
 51:     for (i = 0; i < isize; i++) coors[i] = xmin + hx * (i + istart);
 52:     PetscCall(VecRestoreArray(xcoor, &coors));
 53:   } else if (dim == 2) {
 54:     if (bx == DM_BOUNDARY_PERIODIC) hx = (xmax - xmin) / (M);
 55:     else hx = (xmax - xmin) / (M - 1);
 56:     if (by == DM_BOUNDARY_PERIODIC) hy = (ymax - ymin) / (N);
 57:     else hy = (ymax - ymin) / (N - 1);
 58:     PetscCall(VecGetArray(xcoor, &coors));
 59:     cnt = 0;
 60:     for (j = 0; j < jsize; j++) {
 61:       for (i = 0; i < isize; i++) {
 62:         coors[cnt++] = xmin + hx * (i + istart);
 63:         coors[cnt++] = ymin + hy * (j + jstart);
 64:       }
 65:     }
 66:     PetscCall(VecRestoreArray(xcoor, &coors));
 67:   } else if (dim == 3) {
 68:     if (bx == DM_BOUNDARY_PERIODIC) hx = (xmax - xmin) / (M);
 69:     else hx = (xmax - xmin) / (M - 1);
 70:     if (by == DM_BOUNDARY_PERIODIC) hy = (ymax - ymin) / (N);
 71:     else hy = (ymax - ymin) / (N - 1);
 72:     if (bz == DM_BOUNDARY_PERIODIC) hz_ = (zmax - zmin) / (P);
 73:     else hz_ = (zmax - zmin) / (P - 1);
 74:     PetscCall(VecGetArray(xcoor, &coors));
 75:     cnt = 0;
 76:     for (k = 0; k < ksize; k++) {
 77:       for (j = 0; j < jsize; j++) {
 78:         for (i = 0; i < isize; i++) {
 79:           coors[cnt++] = xmin + hx * (i + istart);
 80:           coors[cnt++] = ymin + hy * (j + jstart);
 81:           coors[cnt++] = zmin + hz_ * (k + kstart);
 82:         }
 83:       }
 84:     }
 85:     PetscCall(VecRestoreArray(xcoor, &coors));
 86:   } else SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Cannot create uniform coordinates for this dimension %" PetscInt_FMT, dim);
 87:   PetscCall(DMSetCoordinates(da, xcoor));
 88:   PetscCall(VecDestroy(&xcoor));
 89:   PetscFunctionReturn(PETSC_SUCCESS);
 90: }

 92: /*
 93:     Allows a user to select a subset of the fields to be drawn by VecView() when the vector comes from a DMDA
 94: */
 95: PetscErrorCode DMDASelectFields(DM da, PetscInt *outfields, PetscInt **fields)
 96: {
 97:   PetscInt  step, ndisplayfields, *displayfields, k, j;
 98:   PetscBool flg;

100:   PetscFunctionBegin;
101:   PetscCall(DMDAGetInfo(da, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &step, NULL, NULL, NULL, NULL, NULL));
102:   PetscCall(PetscMalloc1(step, &displayfields));
103:   for (k = 0; k < step; k++) displayfields[k] = k;
104:   ndisplayfields = step;
105:   PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-draw_fields", displayfields, &ndisplayfields, &flg));
106:   if (!ndisplayfields) ndisplayfields = step;
107:   if (!flg) {
108:     char      **fields;
109:     const char *fieldname;
110:     PetscInt    nfields = step;
111:     PetscCall(PetscMalloc1(step, &fields));
112:     PetscCall(PetscOptionsGetStringArray(NULL, NULL, "-draw_fields_by_name", fields, &nfields, &flg));
113:     if (flg) {
114:       ndisplayfields = 0;
115:       for (k = 0; k < nfields; k++) {
116:         for (j = 0; j < step; j++) {
117:           PetscCall(DMDAGetFieldName(da, j, &fieldname));
118:           PetscCall(PetscStrcmp(fieldname, fields[k], &flg));
119:           if (flg) goto found;
120:         }
121:         SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_USER, "Unknown fieldname %s", fields[k]);
122:       found:
123:         displayfields[ndisplayfields++] = j;
124:       }
125:     }
126:     for (k = 0; k < nfields; k++) PetscCall(PetscFree(fields[k]));
127:     PetscCall(PetscFree(fields));
128:   }
129:   *fields    = displayfields;
130:   *outfields = ndisplayfields;
131:   PetscFunctionReturn(PETSC_SUCCESS);
132: }

134: #include <petscdraw.h>

136: PetscErrorCode VecView_MPI_Draw_DA1d(Vec xin, PetscViewer v)
137: {
138:   DM                  da;
139:   PetscMPIInt         rank, size, tag;
140:   PetscInt            i, n, N, dof, istart, isize, j, nbounds;
141:   MPI_Status          status;
142:   PetscReal           min, max, xmin = 0.0, xmax = 0.0, tmp = 0.0, xgtmp = 0.0;
143:   const PetscScalar  *array, *xg;
144:   PetscDraw           draw;
145:   PetscBool           isnull, useports = PETSC_FALSE, showmarkers = PETSC_FALSE;
146:   MPI_Comm            comm;
147:   PetscDrawAxis       axis;
148:   Vec                 xcoor;
149:   DMBoundaryType      bx;
150:   const char         *tlabel = NULL, *xlabel = NULL;
151:   const PetscReal    *bounds;
152:   PetscInt           *displayfields;
153:   PetscInt            k, ndisplayfields;
154:   PetscBool           hold;
155:   PetscDrawViewPorts *ports = NULL;
156:   PetscViewerFormat   format;

158:   PetscFunctionBegin;
159:   PetscCall(PetscViewerDrawGetDraw(v, 0, &draw));
160:   PetscCall(PetscDrawIsNull(draw, &isnull));
161:   if (isnull) PetscFunctionReturn(PETSC_SUCCESS);
162:   PetscCall(PetscViewerDrawGetBounds(v, &nbounds, &bounds));

164:   PetscCall(VecGetDM(xin, &da));
165:   PetscCheck(da, PetscObjectComm((PetscObject)xin), PETSC_ERR_ARG_WRONG, "Vector not generated from a DMDA");
166:   PetscCall(PetscObjectGetComm((PetscObject)xin, &comm));
167:   PetscCallMPI(MPI_Comm_size(comm, &size));
168:   PetscCallMPI(MPI_Comm_rank(comm, &rank));

170:   PetscCall(PetscOptionsGetBool(NULL, NULL, "-draw_vec_use_markers", &showmarkers, NULL));

172:   PetscCall(DMDAGetInfo(da, NULL, &N, NULL, NULL, NULL, NULL, NULL, &dof, NULL, &bx, NULL, NULL, NULL));
173:   PetscCall(DMDAGetCorners(da, &istart, NULL, NULL, &isize, NULL, NULL));
174:   PetscCall(VecGetArrayRead(xin, &array));
175:   PetscCall(VecGetLocalSize(xin, &n));
176:   n = n / dof;

178:   /* Get coordinates of nodes */
179:   PetscCall(DMGetCoordinates(da, &xcoor));
180:   if (!xcoor) {
181:     PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 0.0, 0.0, 0.0));
182:     PetscCall(DMGetCoordinates(da, &xcoor));
183:   }
184:   PetscCall(VecGetArrayRead(xcoor, &xg));
185:   PetscCall(DMDAGetCoordinateName(da, 0, &xlabel));

187:   /* Determine the min and max coordinate in plot */
188:   if (rank == 0) xmin = PetscRealPart(xg[0]);
189:   if (rank == size - 1) xmax = PetscRealPart(xg[n - 1]);
190:   PetscCallMPI(MPI_Bcast(&xmin, 1, MPIU_REAL, 0, comm));
191:   PetscCallMPI(MPI_Bcast(&xmax, 1, MPIU_REAL, size - 1, comm));

193:   PetscCall(DMDASelectFields(da, &ndisplayfields, &displayfields));
194:   PetscCall(PetscViewerGetFormat(v, &format));
195:   PetscCall(PetscOptionsGetBool(NULL, NULL, "-draw_ports", &useports, NULL));
196:   if (format == PETSC_VIEWER_DRAW_PORTS) useports = PETSC_TRUE;
197:   if (useports) {
198:     PetscCall(PetscViewerDrawGetDraw(v, 0, &draw));
199:     PetscCall(PetscViewerDrawGetDrawAxis(v, 0, &axis));
200:     PetscCall(PetscDrawCheckResizedWindow(draw));
201:     PetscCall(PetscDrawClear(draw));
202:     PetscCall(PetscDrawViewPortsCreate(draw, ndisplayfields, &ports));
203:   }

205:   /* Loop over each field; drawing each in a different window */
206:   for (k = 0; k < ndisplayfields; k++) {
207:     j = displayfields[k];

209:     /* determine the min and max value in plot */
210:     PetscCall(VecStrideMin(xin, j, NULL, &min));
211:     PetscCall(VecStrideMax(xin, j, NULL, &max));
212:     if (j < nbounds) {
213:       min = PetscMin(min, bounds[2 * j]);
214:       max = PetscMax(max, bounds[2 * j + 1]);
215:     }
216:     if (min == max) {
217:       min -= 1.e-5;
218:       max += 1.e-5;
219:     }

221:     if (useports) {
222:       PetscCall(PetscDrawViewPortsSet(ports, k));
223:       PetscCall(DMDAGetFieldName(da, j, &tlabel));
224:     } else {
225:       const char *title;
226:       PetscCall(PetscViewerDrawGetHold(v, &hold));
227:       PetscCall(PetscViewerDrawGetDraw(v, k, &draw));
228:       PetscCall(PetscViewerDrawGetDrawAxis(v, k, &axis));
229:       PetscCall(DMDAGetFieldName(da, j, &title));
230:       if (title) PetscCall(PetscDrawSetTitle(draw, title));
231:       PetscCall(PetscDrawCheckResizedWindow(draw));
232:       if (!hold) PetscCall(PetscDrawClear(draw));
233:     }
234:     PetscCall(PetscDrawAxisSetLabels(axis, tlabel, xlabel, NULL));
235:     PetscCall(PetscDrawAxisSetLimits(axis, xmin, xmax, min, max));
236:     PetscCall(PetscDrawAxisDraw(axis));

238:     /* draw local part of vector */
239:     PetscCall(PetscObjectGetNewTag((PetscObject)xin, &tag));
240:     if (rank < size - 1) { /*send value to right */
241:       PetscCallMPI(MPI_Send((void *)&xg[n - 1], 1, MPIU_REAL, rank + 1, tag, comm));
242:       PetscCallMPI(MPI_Send((void *)&array[j + (n - 1) * dof], 1, MPIU_REAL, rank + 1, tag, comm));
243:     }
244:     if (rank) { /* receive value from left */
245:       PetscCallMPI(MPI_Recv(&xgtmp, 1, MPIU_REAL, rank - 1, tag, comm, &status));
246:       PetscCallMPI(MPI_Recv(&tmp, 1, MPIU_REAL, rank - 1, tag, comm, &status));
247:     }
248:     PetscDrawCollectiveBegin(draw);
249:     if (rank) {
250:       PetscCall(PetscDrawLine(draw, xgtmp, tmp, PetscRealPart(xg[0]), PetscRealPart(array[j]), PETSC_DRAW_RED));
251:       if (showmarkers) PetscCall(PetscDrawPoint(draw, xgtmp, tmp, PETSC_DRAW_BLACK));
252:     }
253:     for (i = 1; i < n; i++) {
254:       PetscCall(PetscDrawLine(draw, PetscRealPart(xg[i - 1]), PetscRealPart(array[j + dof * (i - 1)]), PetscRealPart(xg[i]), PetscRealPart(array[j + dof * i]), PETSC_DRAW_RED));
255:       if (showmarkers) PetscCall(PetscDrawMarker(draw, PetscRealPart(xg[i - 1]), PetscRealPart(array[j + dof * (i - 1)]), PETSC_DRAW_BLACK));
256:     }
257:     if (rank == size - 1) {
258:       if (showmarkers) PetscCall(PetscDrawMarker(draw, PetscRealPart(xg[n - 1]), PetscRealPart(array[j + dof * (n - 1)]), PETSC_DRAW_BLACK));
259:     }
260:     PetscDrawCollectiveEnd(draw);
261:     PetscCall(PetscDrawFlush(draw));
262:     PetscCall(PetscDrawPause(draw));
263:     if (!useports) PetscCall(PetscDrawSave(draw));
264:   }
265:   if (useports) {
266:     PetscCall(PetscViewerDrawGetDraw(v, 0, &draw));
267:     PetscCall(PetscDrawSave(draw));
268:   }

270:   PetscCall(PetscDrawViewPortsDestroy(ports));
271:   PetscCall(PetscFree(displayfields));
272:   PetscCall(VecRestoreArrayRead(xcoor, &xg));
273:   PetscCall(VecRestoreArrayRead(xin, &array));
274:   PetscFunctionReturn(PETSC_SUCCESS);
275: }