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 * (PetscReal)xmax);
24: ya = ymax - (int)(yr * (PetscReal)ymax);
25: xb = (int)(xr * (PetscReal)xmax);
26: yb = ymax - (int)(yl * (PetscReal)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 = (short)XTRANS(draw, XiWin, X1);
176: pt[0].y = (short)YTRANS(draw, XiWin, Y_1);
177: pt[1].x = (short)XTRANS(draw, XiWin, X2);
178: pt[1].y = (short)YTRANS(draw, XiWin, Y2);
179: pt[2].x = (short)XTRANS(draw, XiWin, X3);
180: pt[2].y = (short)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, (int)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] = {0, 0, 0, 0};
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: XMaskEvent(win->disp, ButtonReleaseMask, &report);
469: /* get mouse pointer coordinates */
470: XQueryPointer(win->disp, report.xmotion.window, &root, &child, &root_x, &root_y, &px, &py, &keys_button);
471: /* the user may resize the window before pressing the mouse button */
472: XGetGeometry(win->disp, win->win, &root, &root_x, &root_y, &w, &h, &border, &depth);
473: /* cleanup input event handler and cursor */
474: XSelectInput(win->disp, win->win, NoEventMask);
475: XUndefineCursor(win->disp, win->win);
476: XFreeCursor(win->disp, cursor);
477: XSync(win->disp, False);
479: switch (report.xbutton.button) {
480: case Button1:
481: *button = PETSC_BUTTON_LEFT;
482: break;
483: case Button2:
484: *button = PETSC_BUTTON_CENTER;
485: break;
486: case Button3:
487: *button = PETSC_BUTTON_RIGHT;
488: break;
489: case Button4:
490: *button = PETSC_BUTTON_WHEEL_UP;
491: break;
492: case Button5:
493: *button = PETSC_BUTTON_WHEEL_DOWN;
494: break;
495: }
496: if (report.xbutton.state & ShiftMask) {
497: switch (report.xbutton.button) {
498: case Button1:
499: *button = PETSC_BUTTON_LEFT_SHIFT;
500: break;
501: case Button2:
502: *button = PETSC_BUTTON_CENTER_SHIFT;
503: break;
504: case Button3:
505: *button = PETSC_BUTTON_RIGHT_SHIFT;
506: break;
507: }
508: }
509: xx = ((PetscReal)px) / w;
510: yy = 1 - ((PetscReal)py) / h;
511: if (x_user) *x_user = draw->coor_xl + (xx - draw->port_xl) * (draw->coor_xr - draw->coor_xl) / (draw->port_xr - draw->port_xl);
512: if (y_user) *y_user = draw->coor_yl + (yy - draw->port_yl) * (draw->coor_yr - draw->coor_yl) / (draw->port_yr - draw->port_yl);
513: if (x_phys) *x_phys = xx;
514: if (y_phys) *y_phys = yy;
516: finally:
517: PetscDrawCollectiveEnd(draw);
518: PetscCall(PetscDrawCheckResizedWindow_X(draw));
519: PetscFunctionReturn(PETSC_SUCCESS);
520: }
522: static PetscErrorCode PetscDrawSetVisible_X(PetscDraw draw, PetscBool visible)
523: {
524: PetscDraw_X *Xwin = (PetscDraw_X *)draw->data;
525: PetscMPIInt rank;
527: PetscFunctionBegin;
528: PetscDrawCollectiveBegin(draw);
529: PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)draw), &rank));
530: if (rank == 0) {
531: if (Xwin->win) {
532: if (visible) {
533: XMapWindow(Xwin->disp, Xwin->win);
534: } else {
535: XUnmapWindow(Xwin->disp, Xwin->win);
536: XFlush(Xwin->disp);
537: }
538: }
539: }
540: PetscDrawCollectiveEnd(draw);
541: PetscFunctionReturn(PETSC_SUCCESS);
542: }
544: static PetscErrorCode PetscDrawPause_X(PetscDraw draw)
545: {
546: PetscDraw_X *win = (PetscDraw_X *)draw->data;
548: PetscFunctionBegin;
549: if (!win->win) PetscFunctionReturn(PETSC_SUCCESS);
550: if (draw->pause > 0) PetscCall(PetscSleep(draw->pause));
551: else if (draw->pause == -1) {
552: PetscDrawButton button = PETSC_BUTTON_NONE;
553: PetscCall(PetscDrawGetMouseButton(draw, &button, NULL, NULL, NULL, NULL));
554: if (button == PETSC_BUTTON_CENTER) draw->pause = 0;
555: }
556: PetscFunctionReturn(PETSC_SUCCESS);
557: }
559: static PetscErrorCode PetscDrawDestroy_X(PetscDraw draw)
560: {
561: PetscDraw_X *win = (PetscDraw_X *)draw->data;
563: PetscFunctionBegin;
564: PetscCall(PetscDrawDestroy(&draw->popup));
565: PetscCall(PetscDrawXiClose(win));
566: PetscCall(PetscFree(draw->data));
567: PetscFunctionReturn(PETSC_SUCCESS);
568: }
570: static PetscErrorCode PetscDrawGetSingleton_X(PetscDraw, PetscDraw *);
571: static PetscErrorCode PetscDrawRestoreSingleton_X(PetscDraw, PetscDraw *);
572: PETSC_INTERN PetscErrorCode PetscDrawGetImage_X(PetscDraw, unsigned char[][3], unsigned int *, unsigned int *, unsigned char *[]);
574: 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};
576: static PetscErrorCode PetscDrawGetSingleton_X(PetscDraw draw, PetscDraw *sdraw)
577: {
578: PetscDraw_X *Xwin = (PetscDraw_X *)draw->data, *sXwin;
580: PetscFunctionBegin;
581: PetscCall(PetscDrawCreate(PETSC_COMM_SELF, draw->display, draw->title, draw->x, draw->y, draw->w, draw->h, sdraw));
582: PetscCall(PetscObjectChangeTypeName((PetscObject)*sdraw, PETSC_DRAW_X));
583: (*sdraw)->ops[0] = DvOps;
585: if (draw->popup) PetscCall(PetscDrawGetSingleton(draw->popup, &(*sdraw)->popup));
586: (*sdraw)->pause = draw->pause;
587: (*sdraw)->coor_xl = draw->coor_xl;
588: (*sdraw)->coor_xr = draw->coor_xr;
589: (*sdraw)->coor_yl = draw->coor_yl;
590: (*sdraw)->coor_yr = draw->coor_yr;
591: (*sdraw)->port_xl = draw->port_xl;
592: (*sdraw)->port_xr = draw->port_xr;
593: (*sdraw)->port_yl = draw->port_yl;
594: (*sdraw)->port_yr = draw->port_yr;
596: /* share drawables (windows and/or pixmap) from the parent draw */
597: PetscCall(PetscNew(&sXwin));
598: (*sdraw)->data = (void *)sXwin;
599: PetscCall(PetscDrawXiInit(sXwin, draw->display));
600: if (Xwin->win) {
601: PetscCall(PetscDrawXiQuickWindowFromWindow(sXwin, Xwin->win));
602: sXwin->drw = Xwin->drw; /* XXX If the window is ever resized, this is wrong! */
603: } else if (Xwin->drw) {
604: PetscCall(PetscDrawXiColormap(sXwin));
605: sXwin->drw = Xwin->drw;
606: }
607: PetscCall(PetscDrawXiGetGeometry(sXwin, &sXwin->x, &sXwin->y, &sXwin->w, &sXwin->h));
608: (*sdraw)->x = sXwin->x;
609: (*sdraw)->y = sXwin->y;
610: (*sdraw)->w = sXwin->w;
611: (*sdraw)->h = sXwin->h;
612: PetscFunctionReturn(PETSC_SUCCESS);
613: }
615: static PetscErrorCode PetscDrawRestoreSingleton_X(PetscDraw draw, PetscDraw *sdraw)
616: {
617: PetscFunctionBegin;
618: if (draw->popup && (*sdraw)->popup) {
619: PetscBool isdrawx;
620: PetscDraw_X *pXwin = (PetscDraw_X *)draw->popup->data;
621: PetscDraw_X *sXwin = (PetscDraw_X *)(*sdraw)->popup->data;
622: PetscCall(PetscObjectTypeCompare((PetscObject)draw->popup, PETSC_DRAW_X, &isdrawx));
623: if (!isdrawx) goto finally;
624: PetscCall(PetscObjectTypeCompare((PetscObject)(*sdraw)->popup, PETSC_DRAW_X, &isdrawx));
625: if (!isdrawx) goto finally;
626: if (sXwin->win == pXwin->win) PetscCall(PetscDrawRestoreSingleton(draw->popup, &(*sdraw)->popup));
627: }
628: finally:
629: PetscCall(PetscDrawDestroy(sdraw));
630: PetscFunctionReturn(PETSC_SUCCESS);
631: }
633: static PetscErrorCode PetscDrawXGetDisplaySize_Private(const char name[], int *width, int *height, PetscBool *has_display)
634: {
635: Display *display;
637: PetscFunctionBegin;
638: display = XOpenDisplay(name);
639: if (display) {
640: *has_display = PETSC_TRUE;
641: *width = (int)DisplayWidth(display, DefaultScreen(display));
642: *height = (int)DisplayHeight(display, DefaultScreen(display));
643: XCloseDisplay(display);
644: } else {
645: *has_display = PETSC_FALSE;
646: *width = 0;
647: *height = 0;
648: PetscCall((*PetscErrorPrintf)("Unable to open display on %s\n\
649: Make sure your COMPUTE NODES are authorized to connect\n\
650: to this X server and either your DISPLAY variable\n\
651: is set or you use the -display name option\n",
652: name));
653: }
654: PetscFunctionReturn(PETSC_SUCCESS);
655: }
657: /*MC
658: PETSC_DRAW_X - PETSc graphics device that uses either X windows or its virtual version Xvfb
660: Options Database Keys:
661: + -display <display> - sets the display to use
662: . -x_virtual - forces use of a X virtual display Xvfb that will not display anything but -draw_save will still work.
663: Xvfb is automatically started up in PetscSetDisplay() with this option
664: . -draw_size w,h - percentage of screen (either 1, .5, .3, .25), or size in pixels
665: . -geometry x,y,w,h - set location and size in pixels
666: . -draw_virtual - do not open a window (draw on a pixmap), -draw_save will still work
667: - -draw_double_buffer - avoid window flickering (draw on pixmap and flush to window)
669: Level: beginner
671: .seealso: `PetscDraw`, `PetscDrawOpenX()`, `PetscDrawSetDisplay()`, `PetscDrawSetFromOptions()`
672: M*/
674: PETSC_EXTERN PetscErrorCode PetscDrawCreate_X(PetscDraw draw)
675: {
676: PetscDraw_X *Xwin;
677: PetscMPIInt rank;
678: int x = draw->x, y = draw->y, w = draw->w, h = draw->h;
679: static int xavailable = 0, yavailable = 0, ybottom = 0, xmax = 0, ymax = 0;
680: PetscBool set, dvirtual = PETSC_FALSE, doublebuffer = PETSC_TRUE, has_display;
681: PetscInt xywh[4], osize = 4, nsizes = 2;
682: PetscReal sizes[2] = {.3, .3};
683: static size_t DISPLAY_LENGTH = 265;
685: PetscFunctionBegin;
686: /* get the display variable */
687: if (!draw->display) {
688: PetscCall(PetscMalloc1(DISPLAY_LENGTH, &draw->display));
689: PetscCall(PetscGetDisplay(draw->display, DISPLAY_LENGTH));
690: }
692: /* initialize the display size */
693: if (!xmax) {
694: PetscCall(PetscDrawXGetDisplaySize_Private(draw->display, &xmax, &ymax, &has_display));
695: /* if some processors fail on this and others succeed then this is a problem ! */
696: if (!has_display) {
697: PetscCall((*PetscErrorPrintf)("PETSc unable to use X windows\nproceeding without graphics\n"));
698: PetscCall(PetscDrawSetType(draw, PETSC_DRAW_NULL));
699: PetscFunctionReturn(PETSC_SUCCESS);
700: }
701: }
703: /* allow user to set size of drawable */
704: PetscCall(PetscOptionsGetRealArray(((PetscObject)draw)->options, ((PetscObject)draw)->prefix, "-draw_size", sizes, &nsizes, &set));
705: if (set && nsizes == 1 && sizes[0] > 1.0) sizes[1] = sizes[0];
706: if (set) {
707: if (sizes[0] > 1.0) w = (int)sizes[0];
708: else if (sizes[0] == 1.0) w = PETSC_DRAW_FULL_SIZE;
709: else if (sizes[0] == .5) w = PETSC_DRAW_HALF_SIZE;
710: else if (PetscEqualReal(sizes[0], .3)) w = PETSC_DRAW_THIRD_SIZE; /* sizes == 0.3 will never be true */
711: else if (sizes[0] == .25) w = PETSC_DRAW_QUARTER_SIZE;
712: if (sizes[1] > 1.0) h = (int)sizes[1];
713: else if (sizes[1] == 1.0) h = PETSC_DRAW_FULL_SIZE;
714: else if (sizes[1] == .5) h = PETSC_DRAW_HALF_SIZE;
715: else if (PetscEqualReal(sizes[1], .3)) h = PETSC_DRAW_THIRD_SIZE; /* sizes == 0.3 will never be true */
716: else if (sizes[1] == .25) h = PETSC_DRAW_QUARTER_SIZE;
717: }
718: if (w == PETSC_DECIDE || w == PETSC_DEFAULT) w = draw->w = 300;
719: if (h == PETSC_DECIDE || h == PETSC_DEFAULT) h = draw->h = 300;
720: switch (w) {
721: case PETSC_DRAW_FULL_SIZE:
722: w = draw->w = (xmax - 10);
723: break;
724: case PETSC_DRAW_HALF_SIZE:
725: w = draw->w = (xmax - 20) / 2;
726: break;
727: case PETSC_DRAW_THIRD_SIZE:
728: w = draw->w = (xmax - 30) / 3;
729: break;
730: case PETSC_DRAW_QUARTER_SIZE:
731: w = draw->w = (xmax - 40) / 4;
732: break;
733: }
734: switch (h) {
735: case PETSC_DRAW_FULL_SIZE:
736: h = draw->h = (ymax - 10);
737: break;
738: case PETSC_DRAW_HALF_SIZE:
739: h = draw->h = (ymax - 20) / 2;
740: break;
741: case PETSC_DRAW_THIRD_SIZE:
742: h = draw->h = (ymax - 30) / 3;
743: break;
744: case PETSC_DRAW_QUARTER_SIZE:
745: h = draw->h = (ymax - 40) / 4;
746: break;
747: }
749: PetscCall(PetscOptionsGetBool(((PetscObject)draw)->options, ((PetscObject)draw)->prefix, "-draw_virtual", &dvirtual, NULL));
751: if (!dvirtual) {
752: /* allow user to set location and size of window */
753: xywh[0] = x;
754: xywh[1] = y;
755: xywh[2] = w;
756: xywh[3] = h;
757: PetscCall(PetscOptionsGetIntArray(((PetscObject)draw)->options, ((PetscObject)draw)->prefix, "-geometry", xywh, &osize, NULL));
758: x = (int)xywh[0];
759: y = (int)xywh[1];
760: w = (int)xywh[2];
761: h = (int)xywh[3];
762: if (w == PETSC_DECIDE || w == PETSC_DEFAULT) w = 300;
763: if (h == PETSC_DECIDE || h == PETSC_DEFAULT) h = 300;
764: draw->x = x;
765: draw->y = y;
766: draw->w = w;
767: draw->h = h;
769: if (draw->x == PETSC_DECIDE || draw->y == PETSC_DECIDE) {
770: /*
771: PETSc tries to place windows starting in the upper left corner
772: and moving across to the right.
774: +0,0-------------------------------------------+
775: | Region used so far +xavailable,yavailable |
776: | | |
777: | | |
778: +--------------------- +ybottom |
779: | |
780: | |
781: +----------------------------------------------+xmax,ymax
783: */
784: /* First: can we add it to the right? */
785: if (xavailable + w + 10 <= xmax) {
786: x = xavailable;
787: y = yavailable;
788: ybottom = PetscMax(ybottom, y + h + 30);
789: } else {
790: /* No, so add it below on the left */
791: xavailable = x = 0;
792: yavailable = y = ybottom;
793: ybottom = ybottom + h + 30;
794: }
795: }
796: /* update available region */
797: xavailable = PetscMax(xavailable, x + w + 10);
798: if (xavailable >= xmax) {
799: xavailable = 0;
800: yavailable = yavailable + h + 30;
801: ybottom = yavailable;
802: }
803: if (yavailable >= ymax) {
804: y = 0;
805: yavailable = 0;
806: ybottom = 0;
807: }
809: } /* endif (!dvirtual) */
811: PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)draw), &rank));
812: PetscCheck(rank != 0 || (w > 0 && h > 0), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative window width or height");
814: PetscCall(PetscNew(&Xwin));
815: draw->ops[0] = DvOps;
816: draw->data = (void *)Xwin;
818: PetscCall(PetscDrawXiInit(Xwin, draw->display));
819: if (!dvirtual) {
820: Xwin->x = x;
821: Xwin->y = y;
822: Xwin->w = w;
823: Xwin->h = h;
824: if (rank == 0) PetscCall(PetscDrawXiQuickWindow(Xwin, draw->title, x, y, w, h));
825: PetscCallMPI(MPI_Bcast(&Xwin->win, 1, MPI_UNSIGNED_LONG, 0, PetscObjectComm((PetscObject)draw)));
826: if (rank) PetscCall(PetscDrawXiQuickWindowFromWindow(Xwin, Xwin->win));
827: } else {
828: Xwin->x = 0;
829: Xwin->y = 0;
830: Xwin->w = w;
831: Xwin->h = h;
832: PetscCall(PetscDrawXiColormap(Xwin));
833: if (rank == 0) PetscCall(PetscDrawXiQuickPixmap(Xwin));
834: PetscCallMPI(MPI_Bcast(&Xwin->drw, 1, MPI_UNSIGNED_LONG, 0, PetscObjectComm((PetscObject)draw)));
835: }
836: PetscCall(PetscDrawXiGetGeometry(Xwin, &Xwin->x, &Xwin->y, &Xwin->w, &Xwin->h));
837: draw->x = Xwin->x;
838: draw->y = Xwin->y;
839: draw->w = Xwin->w;
840: draw->h = Xwin->h;
842: PetscCall(PetscOptionsGetBool(((PetscObject)draw)->options, ((PetscObject)draw)->prefix, "-draw_double_buffer", &doublebuffer, NULL));
843: if (doublebuffer) PetscCall(PetscDrawSetDoubleBuffer(draw));
844: PetscFunctionReturn(PETSC_SUCCESS);
845: }