Actual source code: xops.c

  1: /*
  2:     Defines the operations for the X PetscDraw implementation.
  3: */

  5: #include <../src/sys/classes/draw/impls/x/ximpl.h>

  7: /*
  8:      These macros transform from the users coordinates to the  X-window pixel coordinates.
  9: */
 10: #define XTRANS(draw, xwin, x) ((int)(((xwin)->w - 1) * ((draw)->port_xl + (((x - (draw)->coor_xl) * ((draw)->port_xr - (draw)->port_xl)) / ((draw)->coor_xr - (draw)->coor_xl)))))
 11: #define YTRANS(draw, xwin, y) (((xwin)->h - 1) - (int)(((xwin)->h - 1) * ((draw)->port_yl + (((y - (draw)->coor_yl) * ((draw)->port_yr - (draw)->port_yl)) / ((draw)->coor_yr - (draw)->coor_yl)))))

 13: #define ITRANS(draw, xwin, i) ((draw)->coor_xl + (((PetscReal)(i)) * ((draw)->coor_xr - (draw)->coor_xl) / ((xwin)->w - 1) - (draw)->port_xl) / ((draw)->port_xr - (draw)->port_xl))
 14: #define JTRANS(draw, xwin, j) ((draw)->coor_yl + (((PetscReal)(j)) / ((xwin)->h - 1) + (draw)->port_yl - 1) * ((draw)->coor_yr - (draw)->coor_yl) / ((draw)->port_yl - (draw)->port_yr))

 16: static PetscErrorCode PetscDrawSetViewport_X(PetscDraw draw, PetscReal xl, PetscReal yl, PetscReal xr, PetscReal yr)
 17: {
 18:   PetscDraw_X *XiWin = (PetscDraw_X *)draw->data;
 19:   int          xa, ya, xb, yb, xmax = XiWin->w - 1, ymax = XiWin->h - 1;
 20:   XRectangle   box;

 22:   PetscFunctionBegin;
 23:   xa = (int)(xl * xmax);
 24:   ya = ymax - (int)(yr * ymax);
 25:   xb = (int)(xr * xmax);
 26:   yb = ymax - (int)(yl * ymax);
 27:   PetscDrawCollectiveBegin(draw);
 28:   box.x      = (short)xa;
 29:   box.width  = (unsigned short)(xb + 1 - xa);
 30:   box.y      = (short)ya;
 31:   box.height = (unsigned short)(yb + 1 - ya);
 32:   XSetClipRectangles(XiWin->disp, XiWin->gc.set, 0, 0, &box, 1, Unsorted);
 33:   PetscDrawCollectiveEnd(draw);
 34:   PetscFunctionReturn(PETSC_SUCCESS);
 35: }

 37: static PetscErrorCode PetscDrawCoordinateToPixel_X(PetscDraw draw, PetscReal x, PetscReal y, int *i, int *j)
 38: {
 39:   PetscDraw_X *XiWin = (PetscDraw_X *)draw->data;

 41:   PetscFunctionBegin;
 42:   *i = XTRANS(draw, XiWin, x);
 43:   *j = YTRANS(draw, XiWin, y);
 44:   PetscFunctionReturn(PETSC_SUCCESS);
 45: }

 47: static PetscErrorCode PetscDrawPixelToCoordinate_X(PetscDraw draw, int i, int j, PetscReal *x, PetscReal *y)
 48: {
 49:   PetscDraw_X *XiWin = (PetscDraw_X *)draw->data;

 51:   PetscFunctionBegin;
 52:   *x = ITRANS(draw, XiWin, i);
 53:   *y = JTRANS(draw, XiWin, j);
 54:   PetscFunctionReturn(PETSC_SUCCESS);
 55: }

 57: static PetscErrorCode PetscDrawPoint_X(PetscDraw draw, PetscReal x, PetscReal y, int c)
 58: {
 59:   PetscDraw_X *XiWin = (PetscDraw_X *)draw->data;
 60:   int          xx, yy, i, j;

 62:   PetscFunctionBegin;
 63:   xx = XTRANS(draw, XiWin, x);
 64:   yy = YTRANS(draw, XiWin, y);
 65:   PetscDrawXiSetColor(XiWin, c);
 66:   for (i = -1; i < 2; i++) {
 67:     for (j = -1; j < 2; j++) XDrawPoint(XiWin->disp, PetscDrawXiDrawable(XiWin), XiWin->gc.set, xx + i, yy + j);
 68:   }
 69:   PetscFunctionReturn(PETSC_SUCCESS);
 70: }

 72: static PetscErrorCode PetscDrawPointPixel_X(PetscDraw draw, int x, int y, int c)
 73: {
 74:   PetscDraw_X *XiWin = (PetscDraw_X *)draw->data;

 76:   PetscFunctionBegin;
 77:   PetscDrawXiSetColor(XiWin, c);
 78:   XDrawPoint(XiWin->disp, PetscDrawXiDrawable(XiWin), XiWin->gc.set, x, y);
 79:   PetscFunctionReturn(PETSC_SUCCESS);
 80: }

 82: static PetscErrorCode PetscDrawLine_X(PetscDraw draw, PetscReal xl, PetscReal yl, PetscReal xr, PetscReal yr, int cl)
 83: {
 84:   PetscDraw_X *XiWin = (PetscDraw_X *)draw->data;
 85:   int          x_1, y_1, x_2, y_2;

 87:   PetscFunctionBegin;
 88:   PetscDrawXiSetColor(XiWin, cl);
 89:   x_1 = XTRANS(draw, XiWin, xl);
 90:   x_2 = XTRANS(draw, XiWin, xr);
 91:   y_1 = YTRANS(draw, XiWin, yl);
 92:   y_2 = YTRANS(draw, XiWin, yr);
 93:   XDrawLine(XiWin->disp, PetscDrawXiDrawable(XiWin), XiWin->gc.set, x_1, y_1, x_2, y_2);
 94:   PetscFunctionReturn(PETSC_SUCCESS);
 95: }

 97: static PetscErrorCode PetscDrawArrow_X(PetscDraw draw, PetscReal xl, PetscReal yl, PetscReal xr, PetscReal yr, int cl)
 98: {
 99:   PetscDraw_X *XiWin = (PetscDraw_X *)draw->data;
100:   int          x_1, y_1, x_2, y_2;

102:   PetscFunctionBegin;
103:   PetscDrawXiSetColor(XiWin, cl);
104:   x_1 = XTRANS(draw, XiWin, xl);
105:   x_2 = XTRANS(draw, XiWin, xr);
106:   y_1 = YTRANS(draw, XiWin, yl);
107:   y_2 = YTRANS(draw, XiWin, yr);
108:   XDrawLine(XiWin->disp, PetscDrawXiDrawable(XiWin), XiWin->gc.set, x_1, y_1, x_2, y_2);
109:   if (x_1 == x_2 && y_1 == y_2) PetscFunctionReturn(PETSC_SUCCESS);
110:   if (x_1 == x_2 && PetscAbs(y_1 - y_2) > 7) {
111:     if (y_2 > y_1) {
112:       XDrawLine(XiWin->disp, PetscDrawXiDrawable(XiWin), XiWin->gc.set, x_2, y_2, x_2 - 3, y_2 - 3);
113:       XDrawLine(XiWin->disp, PetscDrawXiDrawable(XiWin), XiWin->gc.set, x_2, y_2, x_2 + 3, y_2 - 3);
114:     } else {
115:       XDrawLine(XiWin->disp, PetscDrawXiDrawable(XiWin), XiWin->gc.set, x_2, y_2, x_2 - 3, y_2 + 3);
116:       XDrawLine(XiWin->disp, PetscDrawXiDrawable(XiWin), XiWin->gc.set, x_2, y_2, x_2 + 3, y_2 + 3);
117:     }
118:   }
119:   if (y_1 == y_2 && PetscAbs(x_1 - x_2) > 7) {
120:     if (x_2 > x_1) {
121:       XDrawLine(XiWin->disp, PetscDrawXiDrawable(XiWin), XiWin->gc.set, x_2 - 3, y_2 - 3, x_2, y_2);
122:       XDrawLine(XiWin->disp, PetscDrawXiDrawable(XiWin), XiWin->gc.set, x_2 - 3, y_2 + 3, x_2, y_2);
123:     } else {
124:       XDrawLine(XiWin->disp, PetscDrawXiDrawable(XiWin), XiWin->gc.set, x_2, y_2, x_2 + 3, y_2 - 3);
125:       XDrawLine(XiWin->disp, PetscDrawXiDrawable(XiWin), XiWin->gc.set, x_2, y_2, x_2 + 3, y_2 + 3);
126:     }
127:   }
128:   PetscFunctionReturn(PETSC_SUCCESS);
129: }

131: static PetscErrorCode PetscDrawRectangle_X(PetscDraw draw, PetscReal xl, PetscReal yl, PetscReal xr, PetscReal yr, int c1, int c2, int c3, int c4)
132: {
133:   PetscDraw_X *XiWin = (PetscDraw_X *)draw->data;
134:   int          x, y, w, h, c = (c1 + c2 + c3 + c4) / 4;

136:   PetscFunctionBegin;
137:   PetscDrawXiSetColor(XiWin, c);
138:   x = XTRANS(draw, XiWin, xl);
139:   w = XTRANS(draw, XiWin, xr) + 1 - x;
140:   if (w <= 0) w = 1;
141:   y = YTRANS(draw, XiWin, yr);
142:   h = YTRANS(draw, XiWin, yl) + 1 - y;
143:   if (h <= 0) h = 1;
144:   XFillRectangle(XiWin->disp, PetscDrawXiDrawable(XiWin), XiWin->gc.set, x, y, w, h);
145:   PetscFunctionReturn(PETSC_SUCCESS);
146: }

148: static PetscErrorCode PetscDrawEllipse_X(PetscDraw draw, PetscReal x, PetscReal y, PetscReal a, PetscReal b, int c)
149: {
150:   PetscDraw_X *XiWin = (PetscDraw_X *)draw->data;
151:   int          xA, yA, w, h;

153:   PetscFunctionBegin;
154:   PetscDrawXiSetColor(XiWin, c);
155:   xA = XTRANS(draw, XiWin, x - a / 2);
156:   w  = XTRANS(draw, XiWin, x + a / 2) + 1 - xA;
157:   w  = PetscAbs(w);
158:   yA = YTRANS(draw, XiWin, y + b / 2);
159:   h  = YTRANS(draw, XiWin, y - b / 2) + 1 - yA;
160:   h  = PetscAbs(h);
161:   XFillArc(XiWin->disp, PetscDrawXiDrawable(XiWin), XiWin->gc.set, xA, yA, w, h, 0, 360 * 64);
162:   PetscFunctionReturn(PETSC_SUCCESS);
163: }

165: PETSC_INTERN PetscErrorCode PetscDrawInterpolatedTriangle_X(PetscDraw_X *, int, int, int, int, int, int, int, int, int);

167: static PetscErrorCode PetscDrawTriangle_X(PetscDraw draw, PetscReal X1, PetscReal Y_1, PetscReal X2, PetscReal Y2, PetscReal X3, PetscReal Y3, int c1, int c2, int c3)
168: {
169:   PetscDraw_X *XiWin = (PetscDraw_X *)draw->data;

171:   PetscFunctionBegin;
172:   if (c1 == c2 && c2 == c3) {
173:     XPoint pt[3];
174:     PetscDrawXiSetColor(XiWin, c1);
175:     pt[0].x = XTRANS(draw, XiWin, X1);
176:     pt[0].y = YTRANS(draw, XiWin, Y_1);
177:     pt[1].x = XTRANS(draw, XiWin, X2);
178:     pt[1].y = YTRANS(draw, XiWin, Y2);
179:     pt[2].x = XTRANS(draw, XiWin, X3);
180:     pt[2].y = YTRANS(draw, XiWin, Y3);
181:     XFillPolygon(XiWin->disp, PetscDrawXiDrawable(XiWin), XiWin->gc.set, pt, 3, Convex, CoordModeOrigin);
182:   } else {
183:     int x1, y_1, x2, y2, x3, y3;
184:     x1  = XTRANS(draw, XiWin, X1);
185:     y_1 = YTRANS(draw, XiWin, Y_1);
186:     x2  = XTRANS(draw, XiWin, X2);
187:     y2  = YTRANS(draw, XiWin, Y2);
188:     x3  = XTRANS(draw, XiWin, X3);
189:     y3  = YTRANS(draw, XiWin, Y3);
190:     PetscCall(PetscDrawInterpolatedTriangle_X(XiWin, x1, y_1, c1, x2, y2, c2, x3, y3, c3));
191:   }
192:   PetscFunctionReturn(PETSC_SUCCESS);
193: }

195: static PetscErrorCode PetscDrawStringSetSize_X(PetscDraw draw, PetscReal x, PetscReal y)
196: {
197:   PetscDraw_X *XiWin = (PetscDraw_X *)draw->data;
198:   int          w, h;

200:   PetscFunctionBegin;
201:   w = (int)((XiWin->w) * x * (draw->port_xr - draw->port_xl) / (draw->coor_xr - draw->coor_xl));
202:   h = (int)((XiWin->h) * y * (draw->port_yr - draw->port_yl) / (draw->coor_yr - draw->coor_yl));
203:   PetscCall(PetscFree(XiWin->font));
204:   PetscCall(PetscDrawXiFontFixed(XiWin, w, h, &XiWin->font));
205:   PetscFunctionReturn(PETSC_SUCCESS);
206: }

208: static PetscErrorCode PetscDrawStringGetSize_X(PetscDraw draw, PetscReal *x, PetscReal *y)
209: {
210:   PetscDraw_X *XiWin = (PetscDraw_X *)draw->data;
211:   PetscReal    w, h;

213:   PetscFunctionBegin;
214:   w = XiWin->font->font_w;
215:   h = XiWin->font->font_h;
216:   if (x) *x = w * (draw->coor_xr - draw->coor_xl) / ((XiWin->w) * (draw->port_xr - draw->port_xl));
217:   if (y) *y = h * (draw->coor_yr - draw->coor_yl) / ((XiWin->h) * (draw->port_yr - draw->port_yl));
218:   PetscFunctionReturn(PETSC_SUCCESS);
219: }

221: static PetscErrorCode PetscDrawString_X(PetscDraw draw, PetscReal x, PetscReal y, int c, const char chrs[])
222: {
223:   PetscDraw_X *XiWin = (PetscDraw_X *)draw->data;
224:   int          xx, yy, descent = XiWin->font->font_descent;
225:   size_t       len;
226:   char        *substr;
227:   PetscToken   token;

229:   PetscFunctionBegin;
230:   xx = XTRANS(draw, XiWin, x);
231:   yy = YTRANS(draw, XiWin, y);
232:   PetscDrawXiSetColor(XiWin, c);

234:   PetscCall(PetscTokenCreate(chrs, '\n', &token));
235:   PetscCall(PetscTokenFind(token, &substr));
236:   while (substr) {
237:     PetscCall(PetscStrlen(substr, &len));
238:     XDrawString(XiWin->disp, PetscDrawXiDrawable(XiWin), XiWin->gc.set, xx, yy - descent, substr, len);
239:     yy += XiWin->font->font_h;
240:     PetscCall(PetscTokenFind(token, &substr));
241:   }
242:   PetscCall(PetscTokenDestroy(&token));
243:   PetscFunctionReturn(PETSC_SUCCESS);
244: }

246: static PetscErrorCode PetscDrawStringVertical_X(PetscDraw draw, PetscReal x, PetscReal y, int c, const char text[])
247: {
248:   PetscDraw_X *XiWin = (PetscDraw_X *)draw->data;
249:   int          xx, yy, offset = XiWin->font->font_h - XiWin->font->font_descent;
250:   char         chr[2] = {0, 0};

252:   PetscFunctionBegin;
253:   xx = XTRANS(draw, XiWin, x);
254:   yy = YTRANS(draw, XiWin, y);
255:   PetscDrawXiSetColor(XiWin, c);
256:   while ((chr[0] = *text++)) {
257:     XDrawString(XiWin->disp, PetscDrawXiDrawable(XiWin), XiWin->gc.set, xx, yy + offset, chr, 1);
258:     yy += XiWin->font->font_h;
259:   }
260:   PetscFunctionReturn(PETSC_SUCCESS);
261: }

263: static PetscErrorCode PetscDrawFlush_X(PetscDraw draw)
264: {
265:   PetscDraw_X *XiWin = (PetscDraw_X *)draw->data;
266:   PetscMPIInt  rank;

268:   PetscFunctionBegin;
269:   /* make sure the X server processed requests from all processes */
270:   PetscDrawCollectiveBegin(draw);
271:   XSync(XiWin->disp, False);
272:   PetscDrawCollectiveEnd(draw);
273:   PetscCallMPI(MPI_Barrier(PetscObjectComm((PetscObject)draw)));

275:   /* transfer pixmap contents to window (only the first process does this) */
276:   if (XiWin->drw && XiWin->win) {
277:     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)draw), &rank));
278:     PetscDrawCollectiveBegin(draw);
279:     if (rank == 0) XCopyArea(XiWin->disp, XiWin->drw, XiWin->win, XiWin->gc.set, 0, 0, XiWin->w, XiWin->h, 0, 0);
280:     if (rank == 0) XSync(XiWin->disp, False);
281:     PetscDrawCollectiveEnd(draw);
282:     PetscCallMPI(MPI_Barrier(PetscObjectComm((PetscObject)draw)));
283:   }
284:   PetscFunctionReturn(PETSC_SUCCESS);
285: }

287: static PetscErrorCode PetscDrawClear_X(PetscDraw draw)
288: {
289:   PetscDraw_X *XiWin = (PetscDraw_X *)draw->data;
290:   int          xmax = XiWin->w - 1, ymax = XiWin->h - 1;
291:   PetscReal    xl = draw->port_xl, yl = draw->port_yl;
292:   PetscReal    xr = draw->port_xr, yr = draw->port_yr;
293:   PetscMPIInt  rank;

295:   PetscFunctionBegin;
296:   /* make sure the X server processed requests from all processes */
297:   PetscDrawCollectiveBegin(draw);
298:   XSync(XiWin->disp, False);
299:   PetscDrawCollectiveEnd(draw);
300:   PetscCallMPI(MPI_Barrier(PetscObjectComm((PetscObject)draw)));

302:   /* only the first process handles the clearing business */
303:   PetscDrawCollectiveBegin(draw);
304:   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)draw), &rank));
305:   if (rank == 0) {
306:     int          xa = (int)(xl * xmax), ya = ymax - (int)(yr * ymax);
307:     int          xb = (int)(xr * xmax), yb = ymax - (int)(yl * ymax);
308:     unsigned int w = (unsigned int)(xb + 1 - xa);
309:     unsigned int h = (unsigned int)(yb + 1 - ya);
310:     PetscDrawXiSetPixVal(XiWin, XiWin->background);
311:     XFillRectangle(XiWin->disp, PetscDrawXiDrawable(XiWin), XiWin->gc.set, xa, ya, w, h);
312:     XSync(XiWin->disp, False);
313:   }
314:   PetscDrawCollectiveEnd(draw);
315:   PetscCallMPI(MPI_Barrier(PetscObjectComm((PetscObject)draw)));
316:   PetscFunctionReturn(PETSC_SUCCESS);
317: }

319: static PetscErrorCode PetscDrawSetDoubleBuffer_X(PetscDraw draw)
320: {
321:   PetscDraw_X *win = (PetscDraw_X *)draw->data;
322:   PetscMPIInt  rank;

324:   PetscFunctionBegin;
325:   if (win->drw) PetscFunctionReturn(PETSC_SUCCESS);
326:   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)draw), &rank));

328:   PetscDrawCollectiveBegin(draw);
329:   if (rank == 0) PetscCall(PetscDrawXiQuickPixmap(win));
330:   PetscDrawCollectiveEnd(draw);
331:   PetscCallMPI(MPI_Bcast(&win->drw, 1, MPI_UNSIGNED_LONG, 0, PetscObjectComm((PetscObject)draw)));
332:   PetscFunctionReturn(PETSC_SUCCESS);
333: }

335: static PetscErrorCode PetscDrawGetPopup_X(PetscDraw draw, PetscDraw *popup)
336: {
337:   PetscDraw_X *win = (PetscDraw_X *)draw->data;
338:   PetscBool    flg = PETSC_TRUE;

340:   PetscFunctionBegin;
341:   PetscCall(PetscOptionsGetBool(((PetscObject)draw)->options, ((PetscObject)draw)->prefix, "-draw_popup", &flg, NULL));
342:   if (!flg || !win->win) {
343:     *popup = NULL;
344:     PetscFunctionReturn(PETSC_SUCCESS);
345:   }

347:   PetscCall(PetscDrawCreate(PetscObjectComm((PetscObject)draw), draw->display, NULL, win->x, win->y + win->h + 10, 220, 220, popup));
348:   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*popup, "popup_"));
349:   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)*popup, ((PetscObject)draw)->prefix));
350:   PetscCall(PetscDrawSetType(*popup, PETSC_DRAW_X));
351:   draw->popup = *popup;
352:   PetscFunctionReturn(PETSC_SUCCESS);
353: }

355: static PetscErrorCode PetscDrawSetTitle_X(PetscDraw draw, const char title[])
356: {
357:   PetscDraw_X *win = (PetscDraw_X *)draw->data;
358:   PetscMPIInt  rank;

360:   PetscFunctionBegin;
361:   if (!win->win) PetscFunctionReturn(PETSC_SUCCESS);
362:   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)draw), &rank));
363:   PetscDrawCollectiveBegin(draw);
364:   if (rank == 0) {
365:     size_t        len;
366:     XTextProperty prop;
367:     PetscCall(PetscStrlen(title, &len));
368:     XGetWMName(win->disp, win->win, &prop);
369:     XFree((void *)prop.value);
370:     prop.value  = (unsigned char *)title;
371:     prop.nitems = (long)len;
372:     XSetWMName(win->disp, win->win, &prop);
373:   }
374:   PetscDrawCollectiveEnd(draw);
375:   PetscFunctionReturn(PETSC_SUCCESS);
376: }

378: static PetscErrorCode PetscDrawCheckResizedWindow_X(PetscDraw draw)
379: {
380:   PetscDraw_X *win = (PetscDraw_X *)draw->data;
381:   int          xywh[4];
382:   PetscMPIInt  rank;

384:   PetscFunctionBegin;
385:   if (!win->win) PetscFunctionReturn(PETSC_SUCCESS);
386:   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)draw), &rank));

388:   PetscDrawCollectiveBegin(draw);
389:   if (rank == 0) PetscCall(PetscDrawXiGetGeometry(win, xywh, xywh + 1, xywh + 2, xywh + 3));
390:   PetscDrawCollectiveEnd(draw);
391:   PetscCallMPI(MPI_Bcast(xywh, 4, MPI_INT, 0, PetscObjectComm((PetscObject)draw)));

393:   /* record new window position */
394:   draw->x = win->x = xywh[0];
395:   draw->y = win->y = xywh[1];
396:   if (xywh[2] == win->w && xywh[3] == win->h) PetscFunctionReturn(PETSC_SUCCESS);
397:   /* record new window sizes */
398:   draw->w = win->w = xywh[2];
399:   draw->h = win->h = xywh[3];

401:   /* recreate pixmap (only first processor does this) */
402:   PetscDrawCollectiveBegin(draw);
403:   if (rank == 0 && win->drw) PetscCall(PetscDrawXiQuickPixmap(win));
404:   PetscDrawCollectiveEnd(draw);
405:   PetscCallMPI(MPI_Bcast(&win->drw, 1, MPI_UNSIGNED_LONG, 0, PetscObjectComm((PetscObject)draw)));
406:   /* reset the clipping */
407:   PetscCall(PetscDrawSetViewport_X(draw, draw->port_xl, draw->port_yl, draw->port_xr, draw->port_yr));
408:   PetscFunctionReturn(PETSC_SUCCESS);
409: }

411: static PetscErrorCode PetscDrawResizeWindow_X(PetscDraw draw, int w, int h)
412: {
413:   PetscDraw_X *win = (PetscDraw_X *)draw->data;
414:   PetscMPIInt  rank;

416:   PetscFunctionBegin;
417:   if (w == win->w && h == win->h) PetscFunctionReturn(PETSC_SUCCESS);
418:   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)draw), &rank));

420:   if (win->win) {
421:     PetscDrawCollectiveBegin(draw);
422:     if (rank == 0) PetscCall(PetscDrawXiResizeWindow(win, w, h));
423:     PetscDrawCollectiveEnd(draw);
424:     PetscCall(PetscDrawCheckResizedWindow_X(draw));
425:   } else if (win->drw) {
426:     draw->w = win->w = w;
427:     draw->h = win->h = h;
428:     /* recreate pixmap (only first processor does this) */
429:     PetscDrawCollectiveBegin(draw);
430:     if (rank == 0) PetscCall(PetscDrawXiQuickPixmap(win));
431:     PetscCallMPI(MPI_Bcast(&win->drw, 1, MPI_UNSIGNED_LONG, 0, PetscObjectComm((PetscObject)draw)));
432:     /* reset the clipping */
433:     PetscDrawCollectiveEnd(draw);
434:     PetscCall(PetscDrawSetViewport_X(draw, draw->port_xl, draw->port_yl, draw->port_xr, draw->port_yr));
435:   }
436:   PetscFunctionReturn(PETSC_SUCCESS);
437: }

439: #include <X11/cursorfont.h>

441: static PetscErrorCode PetscDrawGetMouseButton_X(PetscDraw draw, PetscDrawButton *button, PetscReal *x_user, PetscReal *y_user, PetscReal *x_phys, PetscReal *y_phys)
442: {
443:   PetscDraw_X *win = (PetscDraw_X *)draw->data;
444:   Cursor       cursor;
445:   XEvent       report;
446:   Window       root, child;
447:   int          root_x, root_y, px = 0, py = 0;
448:   unsigned int w, h, border, depth;
449:   unsigned int keys_button;
450:   PetscMPIInt  rank;
451:   PetscReal    xx, yy;

453:   PetscFunctionBegin;
454:   *button = PETSC_BUTTON_NONE;
455:   if (!win->win) PetscFunctionReturn(PETSC_SUCCESS);
456:   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)draw), &rank));

458:   PetscDrawCollectiveBegin(draw);
459:   if (rank) goto finally;

461:   /* change cursor to indicate input */
462:   cursor = XCreateFontCursor(win->disp, XC_hand2);
463:   PetscCheck(cursor, PETSC_COMM_SELF, PETSC_ERR_LIB, "Unable to create X cursor");
464:   XDefineCursor(win->disp, win->win, cursor);
465:   /* wait for mouse button events */
466:   XSelectInput(win->disp, win->win, ButtonPressMask | ButtonReleaseMask);
467:   while (XCheckTypedEvent(win->disp, ButtonPress, &report))
468:     ;
469:   XMaskEvent(win->disp, ButtonReleaseMask, &report);
470:   /* get mouse pointer coordinates */
471:   XQueryPointer(win->disp, report.xmotion.window, &root, &child, &root_x, &root_y, &px, &py, &keys_button);
472:   /* the user may resize the window before pressing the mouse button */
473:   XGetGeometry(win->disp, win->win, &root, &root_x, &root_y, &w, &h, &border, &depth);
474:   /* cleanup input event handler and cursor  */
475:   XSelectInput(win->disp, win->win, NoEventMask);
476:   XUndefineCursor(win->disp, win->win);
477:   XFreeCursor(win->disp, cursor);
478:   XSync(win->disp, False);

480:   switch (report.xbutton.button) {
481:   case Button1:
482:     *button = PETSC_BUTTON_LEFT;
483:     break;
484:   case Button2:
485:     *button = PETSC_BUTTON_CENTER;
486:     break;
487:   case Button3:
488:     *button = PETSC_BUTTON_RIGHT;
489:     break;
490:   case Button4:
491:     *button = PETSC_BUTTON_WHEEL_UP;
492:     break;
493:   case Button5:
494:     *button = PETSC_BUTTON_WHEEL_DOWN;
495:     break;
496:   }
497:   if (report.xbutton.state & ShiftMask) {
498:     switch (report.xbutton.button) {
499:     case Button1:
500:       *button = PETSC_BUTTON_LEFT_SHIFT;
501:       break;
502:     case Button2:
503:       *button = PETSC_BUTTON_CENTER_SHIFT;
504:       break;
505:     case Button3:
506:       *button = PETSC_BUTTON_RIGHT_SHIFT;
507:       break;
508:     }
509:   }
510:   xx = ((PetscReal)px) / w;
511:   yy = 1 - ((PetscReal)py) / h;
512:   if (x_user) *x_user = draw->coor_xl + (xx - draw->port_xl) * (draw->coor_xr - draw->coor_xl) / (draw->port_xr - draw->port_xl);
513:   if (y_user) *y_user = draw->coor_yl + (yy - draw->port_yl) * (draw->coor_yr - draw->coor_yl) / (draw->port_yr - draw->port_yl);
514:   if (x_phys) *x_phys = xx;
515:   if (y_phys) *y_phys = yy;

517: finally:
518:   PetscDrawCollectiveEnd(draw);
519:   PetscCall(PetscDrawCheckResizedWindow_X(draw));
520:   PetscFunctionReturn(PETSC_SUCCESS);
521: }

523: static PetscErrorCode PetscDrawSetVisible_X(PetscDraw draw, PetscBool visible)
524: {
525:   PetscDraw_X *Xwin = (PetscDraw_X *)draw->data;
526:   PetscMPIInt  rank;

528:   PetscFunctionBegin;
529:   PetscDrawCollectiveBegin(draw);
530:   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)draw), &rank));
531:   if (rank == 0) {
532:     if (Xwin->win) {
533:       if (visible) {
534:         XMapWindow(Xwin->disp, Xwin->win);
535:       } else {
536:         XUnmapWindow(Xwin->disp, Xwin->win);
537:         XFlush(Xwin->disp);
538:       }
539:     }
540:   }
541:   PetscDrawCollectiveEnd(draw);
542:   PetscFunctionReturn(PETSC_SUCCESS);
543: }

545: static PetscErrorCode PetscDrawPause_X(PetscDraw draw)
546: {
547:   PetscDraw_X *win = (PetscDraw_X *)draw->data;

549:   PetscFunctionBegin;
550:   if (!win->win) PetscFunctionReturn(PETSC_SUCCESS);
551:   if (draw->pause > 0) PetscCall(PetscSleep(draw->pause));
552:   else if (draw->pause == -1) {
553:     PetscDrawButton button = PETSC_BUTTON_NONE;
554:     PetscCall(PetscDrawGetMouseButton(draw, &button, NULL, NULL, NULL, NULL));
555:     if (button == PETSC_BUTTON_CENTER) draw->pause = 0;
556:   }
557:   PetscFunctionReturn(PETSC_SUCCESS);
558: }

560: static PetscErrorCode PetscDrawDestroy_X(PetscDraw draw)
561: {
562:   PetscDraw_X *win = (PetscDraw_X *)draw->data;

564:   PetscFunctionBegin;
565:   PetscCall(PetscDrawDestroy(&draw->popup));
566:   PetscCall(PetscDrawXiClose(win));
567:   PetscCall(PetscFree(draw->data));
568:   PetscFunctionReturn(PETSC_SUCCESS);
569: }

571: static PetscErrorCode       PetscDrawGetSingleton_X(PetscDraw, PetscDraw *);
572: static PetscErrorCode       PetscDrawRestoreSingleton_X(PetscDraw, PetscDraw *);
573: PETSC_INTERN PetscErrorCode PetscDrawGetImage_X(PetscDraw, unsigned char[][3], unsigned int *, unsigned int *, unsigned char *[]);

575: static struct _PetscDrawOps DvOps = {PetscDrawSetDoubleBuffer_X, PetscDrawFlush_X, PetscDrawLine_X, NULL, NULL, PetscDrawPoint_X, NULL, PetscDrawString_X, PetscDrawStringVertical_X, PetscDrawStringSetSize_X, PetscDrawStringGetSize_X, PetscDrawSetViewport_X, PetscDrawClear_X, PetscDrawRectangle_X, PetscDrawTriangle_X, PetscDrawEllipse_X, PetscDrawGetMouseButton_X, PetscDrawPause_X, NULL, NULL, PetscDrawGetPopup_X, PetscDrawSetTitle_X, PetscDrawCheckResizedWindow_X, PetscDrawResizeWindow_X, PetscDrawDestroy_X, NULL, PetscDrawGetSingleton_X, PetscDrawRestoreSingleton_X, NULL, PetscDrawGetImage_X, NULL, PetscDrawArrow_X, PetscDrawCoordinateToPixel_X, PetscDrawPixelToCoordinate_X, PetscDrawPointPixel_X, NULL, PetscDrawSetVisible_X};

577: static PetscErrorCode PetscDrawGetSingleton_X(PetscDraw draw, PetscDraw *sdraw)
578: {
579:   PetscDraw_X *Xwin = (PetscDraw_X *)draw->data, *sXwin;

581:   PetscFunctionBegin;
582:   PetscCall(PetscDrawCreate(PETSC_COMM_SELF, draw->display, draw->title, draw->x, draw->y, draw->w, draw->h, sdraw));
583:   PetscCall(PetscObjectChangeTypeName((PetscObject)*sdraw, PETSC_DRAW_X));
584:   (*sdraw)->ops[0] = DvOps;

586:   if (draw->popup) PetscCall(PetscDrawGetSingleton(draw->popup, &(*sdraw)->popup));
587:   (*sdraw)->pause   = draw->pause;
588:   (*sdraw)->coor_xl = draw->coor_xl;
589:   (*sdraw)->coor_xr = draw->coor_xr;
590:   (*sdraw)->coor_yl = draw->coor_yl;
591:   (*sdraw)->coor_yr = draw->coor_yr;
592:   (*sdraw)->port_xl = draw->port_xl;
593:   (*sdraw)->port_xr = draw->port_xr;
594:   (*sdraw)->port_yl = draw->port_yl;
595:   (*sdraw)->port_yr = draw->port_yr;

597:   /* share drawables (windows and/or pixmap) from the parent draw */
598:   PetscCall(PetscNew(&sXwin));
599:   (*sdraw)->data = (void *)sXwin;
600:   PetscCall(PetscDrawXiInit(sXwin, draw->display));
601:   if (Xwin->win) {
602:     PetscCall(PetscDrawXiQuickWindowFromWindow(sXwin, Xwin->win));
603:     sXwin->drw = Xwin->drw; /* XXX If the window is ever resized, this is wrong! */
604:   } else if (Xwin->drw) {
605:     PetscCall(PetscDrawXiColormap(sXwin));
606:     sXwin->drw = Xwin->drw;
607:   }
608:   PetscCall(PetscDrawXiGetGeometry(sXwin, &sXwin->x, &sXwin->y, &sXwin->w, &sXwin->h));
609:   (*sdraw)->x = sXwin->x;
610:   (*sdraw)->y = sXwin->y;
611:   (*sdraw)->w = sXwin->w;
612:   (*sdraw)->h = sXwin->h;
613:   PetscFunctionReturn(PETSC_SUCCESS);
614: }

616: static PetscErrorCode PetscDrawRestoreSingleton_X(PetscDraw draw, PetscDraw *sdraw)
617: {
618:   PetscFunctionBegin;
619:   if (draw->popup && (*sdraw)->popup) {
620:     PetscBool    isdrawx;
621:     PetscDraw_X *pXwin = (PetscDraw_X *)draw->popup->data;
622:     PetscDraw_X *sXwin = (PetscDraw_X *)(*sdraw)->popup->data;
623:     PetscCall(PetscObjectTypeCompare((PetscObject)draw->popup, PETSC_DRAW_X, &isdrawx));
624:     if (!isdrawx) goto finally;
625:     PetscCall(PetscObjectTypeCompare((PetscObject)(*sdraw)->popup, PETSC_DRAW_X, &isdrawx));
626:     if (!isdrawx) goto finally;
627:     if (sXwin->win == pXwin->win) PetscCall(PetscDrawRestoreSingleton(draw->popup, &(*sdraw)->popup));
628:   }
629: finally:
630:   PetscCall(PetscDrawDestroy(sdraw));
631:   PetscFunctionReturn(PETSC_SUCCESS);
632: }

634: static PetscErrorCode PetscDrawXGetDisplaySize_Private(const char name[], int *width, int *height, PetscBool *has_display)
635: {
636:   Display *display;

638:   PetscFunctionBegin;
639:   display = XOpenDisplay(name);
640:   if (display) {
641:     *has_display = PETSC_TRUE;
642:     *width       = (int)DisplayWidth(display, DefaultScreen(display));
643:     *height      = (int)DisplayHeight(display, DefaultScreen(display));
644:     XCloseDisplay(display);
645:   } else {
646:     *has_display = PETSC_FALSE;
647:     *width       = 0;
648:     *height      = 0;
649:     PetscCall((*PetscErrorPrintf)("Unable to open display on %s\n\
650:     Make sure your COMPUTE NODES are authorized to connect\n\
651:     to this X server and either your DISPLAY variable\n\
652:     is set or you use the -display name option\n",
653:                                   name));
654:   }
655:   PetscFunctionReturn(PETSC_SUCCESS);
656: }

658: /*MC
659:      PETSC_DRAW_X  - PETSc graphics device that uses either X windows or its virtual version Xvfb

661:    Options Database Keys:
662: +  -display <display> - sets the display to use
663: .  -x_virtual - forces use of a X virtual display Xvfb that will not display anything but -draw_save will still work.
664:                 Xvfb is automatically started up in PetscSetDisplay() with this option
665: .  -draw_size w,h - percentage of screen (either 1, .5, .3, .25), or size in pixels
666: .  -geometry x,y,w,h - set location and size in pixels
667: .  -draw_virtual - do not open a window (draw on a pixmap), -draw_save will still work
668: -  -draw_double_buffer - avoid window flickering (draw on pixmap and flush to window)

670:    Level: beginner

672: .seealso: `PetscDraw`, `PetscDrawOpenX()`, `PetscDrawSetDisplay()`, `PetscDrawSetFromOptions()`
673: M*/

675: PETSC_EXTERN PetscErrorCode PetscDrawCreate_X(PetscDraw draw)
676: {
677:   PetscDraw_X  *Xwin;
678:   PetscMPIInt   rank;
679:   int           x = draw->x, y = draw->y, w = draw->w, h = draw->h;
680:   static int    xavailable = 0, yavailable = 0, ybottom = 0, xmax = 0, ymax = 0;
681:   PetscBool     set, dvirtual = PETSC_FALSE, doublebuffer = PETSC_TRUE, has_display;
682:   PetscInt      xywh[4], osize = 4, nsizes = 2;
683:   PetscReal     sizes[2]       = {.3, .3};
684:   static size_t DISPLAY_LENGTH = 265;

686:   PetscFunctionBegin;
687:   /* get the display variable */
688:   if (!draw->display) {
689:     PetscCall(PetscMalloc1(DISPLAY_LENGTH, &draw->display));
690:     PetscCall(PetscGetDisplay(draw->display, DISPLAY_LENGTH));
691:   }

693:   /* initialize the display size */
694:   if (!xmax) {
695:     PetscCall(PetscDrawXGetDisplaySize_Private(draw->display, &xmax, &ymax, &has_display));
696:     /* if some processors fail on this and others succeed then this is a problem ! */
697:     if (!has_display) {
698:       PetscCall((*PetscErrorPrintf)("PETSc unable to use X windows\nproceeding without graphics\n"));
699:       PetscCall(PetscDrawSetType(draw, PETSC_DRAW_NULL));
700:       PetscFunctionReturn(PETSC_SUCCESS);
701:     }
702:   }

704:   /* allow user to set size of drawable */
705:   PetscCall(PetscOptionsGetRealArray(((PetscObject)draw)->options, ((PetscObject)draw)->prefix, "-draw_size", sizes, &nsizes, &set));
706:   if (set && nsizes == 1 && sizes[0] > 1.0) sizes[1] = sizes[0];
707:   if (set) {
708:     if (sizes[0] > 1.0) w = (int)sizes[0];
709:     else if (sizes[0] == 1.0) w = PETSC_DRAW_FULL_SIZE;
710:     else if (sizes[0] == .5) w = PETSC_DRAW_HALF_SIZE;
711:     else if (PetscEqualReal(sizes[0], .3)) w = PETSC_DRAW_THIRD_SIZE; /* sizes == 0.3 will never be true */
712:     else if (sizes[0] == .25) w = PETSC_DRAW_QUARTER_SIZE;
713:     if (sizes[1] > 1.0) h = (int)sizes[1];
714:     else if (sizes[1] == 1.0) h = PETSC_DRAW_FULL_SIZE;
715:     else if (sizes[1] == .5) h = PETSC_DRAW_HALF_SIZE;
716:     else if (PetscEqualReal(sizes[1], .3)) h = PETSC_DRAW_THIRD_SIZE; /* sizes == 0.3 will never be true */
717:     else if (sizes[1] == .25) h = PETSC_DRAW_QUARTER_SIZE;
718:   }
719:   if (w == PETSC_DECIDE || w == PETSC_DEFAULT) w = draw->w = 300;
720:   if (h == PETSC_DECIDE || h == PETSC_DEFAULT) h = draw->h = 300;
721:   switch (w) {
722:   case PETSC_DRAW_FULL_SIZE:
723:     w = draw->w = (xmax - 10);
724:     break;
725:   case PETSC_DRAW_HALF_SIZE:
726:     w = draw->w = (xmax - 20) / 2;
727:     break;
728:   case PETSC_DRAW_THIRD_SIZE:
729:     w = draw->w = (xmax - 30) / 3;
730:     break;
731:   case PETSC_DRAW_QUARTER_SIZE:
732:     w = draw->w = (xmax - 40) / 4;
733:     break;
734:   }
735:   switch (h) {
736:   case PETSC_DRAW_FULL_SIZE:
737:     h = draw->h = (ymax - 10);
738:     break;
739:   case PETSC_DRAW_HALF_SIZE:
740:     h = draw->h = (ymax - 20) / 2;
741:     break;
742:   case PETSC_DRAW_THIRD_SIZE:
743:     h = draw->h = (ymax - 30) / 3;
744:     break;
745:   case PETSC_DRAW_QUARTER_SIZE:
746:     h = draw->h = (ymax - 40) / 4;
747:     break;
748:   }

750:   PetscCall(PetscOptionsGetBool(((PetscObject)draw)->options, ((PetscObject)draw)->prefix, "-draw_virtual", &dvirtual, NULL));

752:   if (!dvirtual) {
753:     /* allow user to set location and size of window */
754:     xywh[0] = x;
755:     xywh[1] = y;
756:     xywh[2] = w;
757:     xywh[3] = h;
758:     PetscCall(PetscOptionsGetIntArray(((PetscObject)draw)->options, ((PetscObject)draw)->prefix, "-geometry", xywh, &osize, NULL));
759:     x = (int)xywh[0];
760:     y = (int)xywh[1];
761:     w = (int)xywh[2];
762:     h = (int)xywh[3];
763:     if (w == PETSC_DECIDE || w == PETSC_DEFAULT) w = 300;
764:     if (h == PETSC_DECIDE || h == PETSC_DEFAULT) h = 300;
765:     draw->x = x;
766:     draw->y = y;
767:     draw->w = w;
768:     draw->h = h;

770:     if (draw->x == PETSC_DECIDE || draw->y == PETSC_DECIDE) {
771:       /*
772:        PETSc tries to place windows starting in the upper left corner
773:         and moving across to the right.

775:        +0,0-------------------------------------------+
776:        |  Region used so far  +xavailable,yavailable  |
777:        |                      |                       |
778:        |                      |                       |
779:        +--------------------- +ybottom                |
780:        |                                              |
781:        |                                              |
782:        +----------------------------------------------+xmax,ymax

784:       */
785:       /*  First: can we add it to the right? */
786:       if (xavailable + w + 10 <= xmax) {
787:         x       = xavailable;
788:         y       = yavailable;
789:         ybottom = PetscMax(ybottom, y + h + 30);
790:       } else {
791:         /* No, so add it below on the left */
792:         xavailable = x = 0;
793:         yavailable = y = ybottom;
794:         ybottom        = ybottom + h + 30;
795:       }
796:     }
797:     /* update available region */
798:     xavailable = PetscMax(xavailable, x + w + 10);
799:     if (xavailable >= xmax) {
800:       xavailable = 0;
801:       yavailable = yavailable + h + 30;
802:       ybottom    = yavailable;
803:     }
804:     if (yavailable >= ymax) {
805:       y          = 0;
806:       yavailable = 0;
807:       ybottom    = 0;
808:     }

810:   } /* endif (!dvirtual) */

812:   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)draw), &rank));
813:   PetscCheck(rank != 0 || (w > 0 && h > 0), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative window width or height");

815:   PetscCall(PetscNew(&Xwin));
816:   draw->ops[0] = DvOps;
817:   draw->data   = (void *)Xwin;

819:   PetscCall(PetscDrawXiInit(Xwin, draw->display));
820:   if (!dvirtual) {
821:     Xwin->x = x;
822:     Xwin->y = y;
823:     Xwin->w = w;
824:     Xwin->h = h;
825:     if (rank == 0) PetscCall(PetscDrawXiQuickWindow(Xwin, draw->title, x, y, w, h));
826:     PetscCallMPI(MPI_Bcast(&Xwin->win, 1, MPI_UNSIGNED_LONG, 0, PetscObjectComm((PetscObject)draw)));
827:     if (rank) PetscCall(PetscDrawXiQuickWindowFromWindow(Xwin, Xwin->win));
828:   } else {
829:     Xwin->x = 0;
830:     Xwin->y = 0;
831:     Xwin->w = w;
832:     Xwin->h = h;
833:     PetscCall(PetscDrawXiColormap(Xwin));
834:     if (rank == 0) PetscCall(PetscDrawXiQuickPixmap(Xwin));
835:     PetscCallMPI(MPI_Bcast(&Xwin->drw, 1, MPI_UNSIGNED_LONG, 0, PetscObjectComm((PetscObject)draw)));
836:   }
837:   PetscCall(PetscDrawXiGetGeometry(Xwin, &Xwin->x, &Xwin->y, &Xwin->w, &Xwin->h));
838:   draw->x = Xwin->x;
839:   draw->y = Xwin->y;
840:   draw->w = Xwin->w;
841:   draw->h = Xwin->h;

843:   PetscCall(PetscOptionsGetBool(((PetscObject)draw)->options, ((PetscObject)draw)->prefix, "-draw_double_buffer", &doublebuffer, NULL));
844:   if (doublebuffer) PetscCall(PetscDrawSetDoubleBuffer(draw));
845:   PetscFunctionReturn(PETSC_SUCCESS);
846: }