Actual source code: dviewp.c

  1: #include <petsc/private/drawimpl.h>

  3: /*@
  4:   PetscDrawSetViewPort - Sets the portion of the window (page) to which draw
  5:   routines will write.

  7:   Collective

  9:   Input Parameters:
 10: + xl   - the horizontal coordinate of the lower left corner of the subwindow.
 11: . yl   - the vertical coordinate of the lower left corner of the subwindow.
 12: . xr   - the horizontal coordinate of the upper right corner of the subwindow.
 13: . yr   - the vertical coordinate of the upper right corner of the subwindow.
 14: - draw - the drawing context

 16:   Level: advanced

 18:   Notes:
 19:   These numbers must always be between 0.0 and 1.0.

 21:   Lower left corner is (0,0).

 23: .seealso: `PetscDrawGetViewPort()`, `PetscDraw`, `PetscDrawSplitViewPort()`, `PetscDrawViewPortsCreate()`
 24: @*/
 25: PetscErrorCode PetscDrawSetViewPort(PetscDraw draw, PetscReal xl, PetscReal yl, PetscReal xr, PetscReal yr)
 26: {
 27:   PetscFunctionBegin;
 29:   PetscCheck(xl >= 0.0 && xr <= 1.0 && yl >= 0.0 && yr <= 1.0 && xr > xl && yr > yl, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "ViewPort values must be >= 0 and <= 1: Instead %g %g %g %g", (double)xl, (double)yl, (double)xr, (double)yr);
 30:   draw->port_xl = xl;
 31:   draw->port_yl = yl;
 32:   draw->port_xr = xr;
 33:   draw->port_yr = yr;
 34:   PetscTryTypeMethod(draw, setviewport, xl, yl, xr, yr);
 35:   PetscFunctionReturn(PETSC_SUCCESS);
 36: }

 38: /*@
 39:   PetscDrawGetViewPort - Gets the portion of the window (page) to which draw
 40:   routines will write.

 42:   Collective

 44:   Input Parameter:
 45: . draw - the drawing context

 47:   Output Parameters:
 48: + xl - the horizontal coordinate of the lower left corner of the subwindow.
 49: . yl - the vertical coordinate of the lower left corner of the subwindow.
 50: . xr - the horizontal coordinate of the upper right corner of the subwindow.
 51: - yr - the vertical coordinate of the upper right corner of the subwindow.

 53:   Level: advanced

 55:   Notes:
 56:   These numbers must always be between 0.0 and 1.0.

 58:   Lower left corner is (0,0).

 60: .seealso: `PetscDraw`, `PetscDrawSplitViewPort()`, `PetscDrawSetViewPort()`
 61: @*/
 62: PetscErrorCode PetscDrawGetViewPort(PetscDraw draw, PetscReal *xl, PetscReal *yl, PetscReal *xr, PetscReal *yr)
 63: {
 64:   PetscFunctionBegin;
 66:   PetscAssertPointer(xl, 2);
 67:   PetscAssertPointer(yl, 3);
 68:   PetscAssertPointer(xr, 4);
 69:   PetscAssertPointer(yr, 5);
 70:   *xl = draw->port_xl;
 71:   *yl = draw->port_yl;
 72:   *xr = draw->port_xr;
 73:   *yr = draw->port_yr;
 74:   PetscFunctionReturn(PETSC_SUCCESS);
 75: }

 77: /*@
 78:   PetscDrawSplitViewPort - Splits a window shared by several processes into smaller
 79:   view ports. One for each process.

 81:   Collective

 83:   Input Parameter:
 84: . draw - the drawing context

 86:   Level: advanced

 88: .seealso: `PetscDrawDivideViewPort()`, `PetscDrawSetViewPort()`
 89: @*/
 90: PetscErrorCode PetscDrawSplitViewPort(PetscDraw draw)
 91: {
 92:   PetscMPIInt rank, size;
 93:   PetscInt    n;
 94:   PetscBool   isnull;
 95:   PetscReal   xl, xr, yl, yr, h;

 97:   PetscFunctionBegin;
 99:   PetscCall(PetscDrawIsNull(draw, &isnull));
100:   if (isnull) PetscFunctionReturn(PETSC_SUCCESS);
101:   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)draw), &rank));
102:   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)draw), &size));

104:   n = (PetscInt)(.1 + PetscSqrtReal((PetscReal)size));
105:   while (n * n < size) n++;

107:   h  = 1.0 / n;
108:   xl = (rank % n) * h;
109:   xr = xl + h;
110:   yl = (rank / n) * h;
111:   yr = yl + h;

113:   PetscDrawCollectiveBegin(draw);
114:   PetscCall(PetscDrawLine(draw, xl, yl, xl, yr, PETSC_DRAW_BLACK));
115:   PetscCall(PetscDrawLine(draw, xl, yr, xr, yr, PETSC_DRAW_BLACK));
116:   PetscCall(PetscDrawLine(draw, xr, yr, xr, yl, PETSC_DRAW_BLACK));
117:   PetscCall(PetscDrawLine(draw, xr, yl, xl, yl, PETSC_DRAW_BLACK));
118:   PetscDrawCollectiveEnd(draw);
119:   PetscCall(PetscDrawFlush(draw));

121:   draw->port_xl = xl + .05 * h;
122:   draw->port_xr = xr - .05 * h;
123:   draw->port_yl = yl + .05 * h;
124:   draw->port_yr = yr - .05 * h;

126:   PetscTryTypeMethod(draw, setviewport, xl, yl, xr, yr);
127:   PetscFunctionReturn(PETSC_SUCCESS);
128: }

130: /*@C
131:   PetscDrawViewPortsCreate - Splits a window into smaller view ports. Each processor shares all the viewports.

133:   Collective

135:   Input Parameters:
136: + draw   - the drawing context
137: - nports - the number of ports

139:   Output Parameter:
140: . newports - a `PetscDrawViewPorts` context (C structure)

142:   Options Database Key:
143: . -draw_ports - display multiple fields in the same window with PetscDrawPorts() instead of in separate windows

145:   Level: advanced

147:   Fortran Note:
148:   No Fortran support since `PetscDrawViewPorts` is a C struct

150: .seealso: `PetscDrawSplitViewPort()`, `PetscDrawSetViewPort()`, `PetscDrawViewPortsSet()`, `PetscDrawViewPortsDestroy()`
151: @*/
152: PetscErrorCode PetscDrawViewPortsCreate(PetscDraw draw, PetscInt nports, PetscDrawViewPorts **newports)
153: {
154:   PetscDrawViewPorts *ports;
155:   PetscInt            i, n;
156:   PetscBool           isnull;
157:   PetscMPIInt         rank;
158:   PetscReal          *xl, *xr, *yl, *yr, h;

160:   PetscFunctionBegin;
162:   PetscCheck(nports >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of divisions must be positive: %" PetscInt_FMT, nports);
163:   PetscAssertPointer(newports, 3);
164:   PetscCall(PetscDrawIsNull(draw, &isnull));
165:   if (isnull) {
166:     *newports = NULL;
167:     PetscFunctionReturn(PETSC_SUCCESS);
168:   }
169:   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)draw), &rank));

171:   PetscCall(PetscNew(&ports));
172:   *newports     = ports;
173:   ports->draw   = draw;
174:   ports->nports = nports;
175:   PetscCall(PetscObjectReference((PetscObject)draw));
176:   /* save previous drawport of window */
177:   PetscCall(PetscDrawGetViewPort(draw, &ports->port_xl, &ports->port_yl, &ports->port_xr, &ports->port_yr));

179:   n = (PetscInt)(.1 + PetscSqrtReal((PetscReal)nports));
180:   while (n * n < nports) n++;
181:   h = 1.0 / n;

183:   PetscCall(PetscMalloc4(n * n, &xl, n * n, &xr, n * n, &yl, n * n, &yr));
184:   ports->xl = xl;
185:   ports->xr = xr;
186:   ports->yl = yl;
187:   ports->yr = yr;

189:   PetscCall(PetscDrawSetCoordinates(draw, 0.0, 0.0, 1.0, 1.0));
190:   PetscDrawCollectiveBegin(draw);
191:   for (i = 0; i < n * n; i++) {
192:     xl[i] = (i % n) * h;
193:     xr[i] = xl[i] + h;
194:     yl[i] = (i / n) * h;
195:     yr[i] = yl[i] + h;

197:     if (rank == 0) {
198:       PetscCall(PetscDrawLine(draw, xl[i], yl[i], xl[i], yr[i], PETSC_DRAW_BLACK));
199:       PetscCall(PetscDrawLine(draw, xl[i], yr[i], xr[i], yr[i], PETSC_DRAW_BLACK));
200:       PetscCall(PetscDrawLine(draw, xr[i], yr[i], xr[i], yl[i], PETSC_DRAW_BLACK));
201:       PetscCall(PetscDrawLine(draw, xr[i], yl[i], xl[i], yl[i], PETSC_DRAW_BLACK));
202:     }

204:     xl[i] += .05 * h;
205:     xr[i] -= .05 * h;
206:     yl[i] += .05 * h;
207:     yr[i] -= .05 * h;
208:   }
209:   PetscDrawCollectiveEnd(draw);
210:   PetscCall(PetscDrawFlush(draw));
211:   PetscFunctionReturn(PETSC_SUCCESS);
212: }

214: /*@C
215:   PetscDrawViewPortsCreateRect - Splits a window into smaller
216:   view ports. Each processor shares all the viewports. The number
217:   of views in the x- and y-directions is specified.

219:   Collective

221:   Input Parameters:
222: + draw - the drawing context
223: . nx   - the number of x divisions
224: - ny   - the number of y divisions

226:   Output Parameter:
227: . newports - a `PetscDrawViewPorts` context (C structure)

229:   Level: advanced

231:   Fortran Note:
232:   No Fortran support since `PetscDrawViewPorts` is a C struct

234: .seealso: `PetscDrawSplitViewPort()`, `PetscDrawSetViewPort()`, `PetscDrawViewPortsSet()`, `PetscDrawViewPortsDestroy()`, `PetscDrawViewPorts`
235: @*/
236: PetscErrorCode PetscDrawViewPortsCreateRect(PetscDraw draw, PetscInt nx, PetscInt ny, PetscDrawViewPorts **newports)
237: {
238:   PetscDrawViewPorts *ports;
239:   PetscReal          *xl, *xr, *yl, *yr, hx, hy;
240:   PetscInt            i, j, k, n;
241:   PetscBool           isnull;
242:   PetscMPIInt         rank;

244:   PetscFunctionBegin;
246:   PetscCheck(nx >= 1 && ny >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of divisions must be positive: %" PetscInt_FMT " x %" PetscInt_FMT, nx, ny);
247:   PetscAssertPointer(newports, 4);
248:   PetscCall(PetscDrawIsNull(draw, &isnull));
249:   if (isnull) {
250:     *newports = NULL;
251:     PetscFunctionReturn(PETSC_SUCCESS);
252:   }
253:   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)draw), &rank));

255:   n  = nx * ny;
256:   hx = 1.0 / nx;
257:   hy = 1.0 / ny;
258:   PetscCall(PetscNew(&ports));
259:   *newports     = ports;
260:   ports->draw   = draw;
261:   ports->nports = n;
262:   PetscCall(PetscObjectReference((PetscObject)draw));
263:   /* save previous drawport of window */
264:   PetscCall(PetscDrawGetViewPort(draw, &ports->port_xl, &ports->port_yl, &ports->port_xr, &ports->port_yr));

266:   PetscCall(PetscMalloc4(n, &xl, n, &xr, n, &yl, n, &yr));
267:   ports->xr = xr;
268:   ports->xl = xl;
269:   ports->yl = yl;
270:   ports->yr = yr;

272:   PetscCall(PetscDrawSetCoordinates(draw, 0.0, 0.0, 1.0, 1.0));
273:   PetscDrawCollectiveBegin(draw);
274:   for (i = 0; i < nx; i++) {
275:     for (j = 0; j < ny; j++) {
276:       k = j * nx + i;

278:       xl[k] = i * hx;
279:       xr[k] = xl[k] + hx;
280:       yl[k] = j * hy;
281:       yr[k] = yl[k] + hy;

283:       if (rank == 0) {
284:         PetscCall(PetscDrawLine(draw, xl[k], yl[k], xl[k], yr[k], PETSC_DRAW_BLACK));
285:         PetscCall(PetscDrawLine(draw, xl[k], yr[k], xr[k], yr[k], PETSC_DRAW_BLACK));
286:         PetscCall(PetscDrawLine(draw, xr[k], yr[k], xr[k], yl[k], PETSC_DRAW_BLACK));
287:         PetscCall(PetscDrawLine(draw, xr[k], yl[k], xl[k], yl[k], PETSC_DRAW_BLACK));
288:       }

290:       xl[k] += .05 * hx;
291:       xr[k] -= .05 * hx;
292:       yl[k] += .05 * hy;
293:       yr[k] -= .05 * hy;
294:     }
295:   }
296:   PetscDrawCollectiveEnd(draw);
297:   PetscCall(PetscDrawFlush(draw));
298:   PetscFunctionReturn(PETSC_SUCCESS);
299: }

301: /*@C
302:   PetscDrawViewPortsDestroy - frees a `PetscDrawViewPorts` object

304:   Collective on the `PetscDraw` inside `ports`

306:   Input Parameter:
307: . ports - the `PetscDrawViewPorts` object

309:   Level: advanced

311:   Fortran Note:
312:   No Fortran support since `PetscDrawViewPorts` is a C struct

314: .seealso: `PetscDrawViewPorts`, `PetscDrawSplitViewPort()`, `PetscDrawSetViewPort()`, `PetscDrawViewPortsSet()`, `PetscDrawViewPortsCreate()`
315: @*/
316: PetscErrorCode PetscDrawViewPortsDestroy(PetscDrawViewPorts *ports)
317: {
318:   PetscFunctionBegin;
319:   if (!ports) PetscFunctionReturn(PETSC_SUCCESS);
320:   PetscAssertPointer(ports, 1);
321:   /* reset Drawport of Window back to previous value */
322:   PetscCall(PetscDrawSetViewPort(ports->draw, ports->port_xl, ports->port_yl, ports->port_xr, ports->port_yr));
323:   PetscCall(PetscDrawDestroy(&ports->draw));
324:   PetscCall(PetscFree4(ports->xl, ports->xr, ports->yl, ports->yr));
325:   PetscCall(PetscFree(ports));
326:   PetscFunctionReturn(PETSC_SUCCESS);
327: }

329: /*@C
330:   PetscDrawViewPortsSet - sets a draw object to use a particular subport

332:   Logically Collective on the `PetscDraw` inside `ports`

334:   Input Parameters:
335: + ports - the `PetscDrawViewPorts` object
336: - port  - the port number, from 0 to nports-1

338:   Level: advanced

340:   Fortran Note:
341:   No Fortran support since `PetscDrawViewPorts` is a C struct

343: .seealso: `PetscDrawViewPorts`, `PetscDrawSplitViewPort()`, `PetscDrawSetViewPort()`, `PetscDrawViewPortsDestroy()`, `PetscDrawViewPortsCreate()`
344: @*/
345: PetscErrorCode PetscDrawViewPortsSet(PetscDrawViewPorts *ports, PetscInt port)
346: {
347:   PetscFunctionBegin;
348:   if (!ports) PetscFunctionReturn(PETSC_SUCCESS);
349:   PetscAssertPointer(ports, 1);
350:   PetscCheck(port >= 0 && (port <= ports->nports - 1), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Port is out of range requested %" PetscInt_FMT " from 0 to %" PetscInt_FMT, port, ports->nports - 1);
351:   PetscCall(PetscDrawSetViewPort(ports->draw, ports->xl[port], ports->yl[port], ports->xr[port], ports->yr[port]));
352:   PetscFunctionReturn(PETSC_SUCCESS);
353: }