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: xa = (int)(xl*xmax); ya = ymax - (int)(yr*ymax);
23: xb = (int)(xr*xmax); yb = ymax - (int)(yl*ymax);
24: PetscDrawCollectiveBegin(draw);
25: box.x = (short)xa; box.width = (unsigned short)(xb + 1 - xa);
26: box.y = (short)ya; box.height = (unsigned short)(yb + 1 - ya);
27: XSetClipRectangles(XiWin->disp,XiWin->gc.set,0,0,&box,1,Unsorted);
28: PetscDrawCollectiveEnd(draw);
29: return 0;
30: }
32: static PetscErrorCode PetscDrawCoordinateToPixel_X(PetscDraw draw,PetscReal x,PetscReal y,int *i,int *j)
33: {
34: PetscDraw_X *XiWin = (PetscDraw_X*)draw->data;
36: *i = XTRANS(draw,XiWin,x);
37: *j = YTRANS(draw,XiWin,y);
38: return 0;
39: }
41: static PetscErrorCode PetscDrawPixelToCoordinate_X(PetscDraw draw,int i,int j,PetscReal *x,PetscReal *y)
42: {
43: PetscDraw_X *XiWin = (PetscDraw_X*)draw->data;
45: *x = ITRANS(draw,XiWin,i);
46: *y = JTRANS(draw,XiWin,j);
47: return 0;
48: }
50: static PetscErrorCode PetscDrawPoint_X(PetscDraw draw,PetscReal x,PetscReal y,int c)
51: {
52: PetscDraw_X *XiWin = (PetscDraw_X*)draw->data;
53: int xx,yy,i,j;
55: xx = XTRANS(draw,XiWin,x);
56: yy = YTRANS(draw,XiWin,y);
57: PetscDrawXiSetColor(XiWin,c);
58: for (i=-1; i<2; i++) {
59: for (j=-1; j<2; j++) {
60: XDrawPoint(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,xx+i,yy+j);
61: }
62: }
63: return 0;
64: }
66: static PetscErrorCode PetscDrawPointPixel_X(PetscDraw draw,int x,int y,int c)
67: {
68: PetscDraw_X *XiWin = (PetscDraw_X*)draw->data;
70: PetscDrawXiSetColor(XiWin,c);
71: XDrawPoint(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x,y);
72: return 0;
73: }
75: static PetscErrorCode PetscDrawLine_X(PetscDraw draw,PetscReal xl,PetscReal yl,PetscReal xr,PetscReal yr,int cl)
76: {
77: PetscDraw_X *XiWin = (PetscDraw_X*)draw->data;
78: int x_1,y_1,x_2,y_2;
80: PetscDrawXiSetColor(XiWin,cl);
81: x_1 = XTRANS(draw,XiWin,xl); x_2 = XTRANS(draw,XiWin,xr);
82: y_1 = YTRANS(draw,XiWin,yl); y_2 = YTRANS(draw,XiWin,yr);
83: XDrawLine(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x_1,y_1,x_2,y_2);
84: return 0;
85: }
87: static PetscErrorCode PetscDrawArrow_X(PetscDraw draw,PetscReal xl,PetscReal yl,PetscReal xr,PetscReal yr,int cl)
88: {
89: PetscDraw_X *XiWin = (PetscDraw_X*)draw->data;
90: int x_1,y_1,x_2,y_2;
92: PetscDrawXiSetColor(XiWin,cl);
93: x_1 = XTRANS(draw,XiWin,xl); x_2 = XTRANS(draw,XiWin,xr);
94: y_1 = YTRANS(draw,XiWin,yl); y_2 = YTRANS(draw,XiWin,yr);
95: XDrawLine(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x_1,y_1,x_2,y_2);
96: if (x_1 == x_2 && y_1 == y_2) return 0;
97: if (x_1 == x_2 && PetscAbs(y_1 - y_2) > 7) {
98: if (y_2 > y_1) {
99: XDrawLine(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x_2,y_2,x_2-3,y_2-3);
100: XDrawLine(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x_2,y_2,x_2+3,y_2-3);
101: } else {
102: XDrawLine(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x_2,y_2,x_2-3,y_2+3);
103: XDrawLine(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x_2,y_2,x_2+3,y_2+3);
104: }
105: }
106: if (y_1 == y_2 && PetscAbs(x_1 - x_2) > 7) {
107: if (x_2 > x_1) {
108: XDrawLine(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x_2-3,y_2-3,x_2,y_2);
109: XDrawLine(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x_2-3,y_2+3,x_2,y_2);
110: } else {
111: XDrawLine(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x_2,y_2,x_2+3,y_2-3);
112: XDrawLine(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x_2,y_2,x_2+3,y_2+3);
113: }
114: }
115: return 0;
116: }
118: static PetscErrorCode PetscDrawRectangle_X(PetscDraw draw,PetscReal xl,PetscReal yl,PetscReal xr,PetscReal yr,int c1,int c2,int c3,int c4)
119: {
120: PetscDraw_X *XiWin = (PetscDraw_X*)draw->data;
121: int x,y,w,h,c = (c1 + c2 + c3 + c4)/4;
123: PetscDrawXiSetColor(XiWin,c);
124: x = XTRANS(draw,XiWin,xl); w = XTRANS(draw,XiWin,xr) + 1 - x; if (w <= 0) w = 1;
125: y = YTRANS(draw,XiWin,yr); h = YTRANS(draw,XiWin,yl) + 1 - y; if (h <= 0) h = 1;
126: XFillRectangle(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x,y,w,h);
127: return 0;
128: }
130: static PetscErrorCode PetscDrawEllipse_X(PetscDraw draw,PetscReal x,PetscReal y,PetscReal a,PetscReal b,int c)
131: {
132: PetscDraw_X *XiWin = (PetscDraw_X*)draw->data;
133: int xA,yA,w,h;
135: PetscDrawXiSetColor(XiWin, c);
136: xA = XTRANS(draw,XiWin, x - a/2); w = XTRANS(draw,XiWin, x + a/2) + 1 - xA; w = PetscAbs(w);
137: yA = YTRANS(draw,XiWin, y + b/2); h = YTRANS(draw,XiWin, y - b/2) + 1 - yA; h = PetscAbs(h);
138: XFillArc(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,xA,yA,w,h,0,360*64);
139: return 0;
140: }
142: PETSC_INTERN PetscErrorCode PetscDrawInterpolatedTriangle_X(PetscDraw_X*,int,int,int,int,int,int,int,int,int);
144: 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)
145: {
146: PetscDraw_X *XiWin = (PetscDraw_X*)draw->data;
148: if (c1 == c2 && c2 == c3) {
149: XPoint pt[3];
150: PetscDrawXiSetColor(XiWin,c1);
151: pt[0].x = XTRANS(draw,XiWin,X1);
152: pt[0].y = YTRANS(draw,XiWin,Y_1);
153: pt[1].x = XTRANS(draw,XiWin,X2);
154: pt[1].y = YTRANS(draw,XiWin,Y2);
155: pt[2].x = XTRANS(draw,XiWin,X3);
156: pt[2].y = YTRANS(draw,XiWin,Y3);
157: XFillPolygon(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,pt,3,Convex,CoordModeOrigin);
158: } else {
159: int x1,y_1,x2,y2,x3,y3;
160: x1 = XTRANS(draw,XiWin,X1);
161: y_1 = YTRANS(draw,XiWin,Y_1);
162: x2 = XTRANS(draw,XiWin,X2);
163: y2 = YTRANS(draw,XiWin,Y2);
164: x3 = XTRANS(draw,XiWin,X3);
165: y3 = YTRANS(draw,XiWin,Y3);
166: PetscDrawInterpolatedTriangle_X(XiWin,x1,y_1,c1,x2,y2,c2,x3,y3,c3);
167: }
168: return 0;
169: }
171: static PetscErrorCode PetscDrawStringSetSize_X(PetscDraw draw,PetscReal x,PetscReal y)
172: {
173: PetscDraw_X *XiWin = (PetscDraw_X*)draw->data;
174: int w,h;
176: w = (int)((XiWin->w)*x*(draw->port_xr - draw->port_xl)/(draw->coor_xr - draw->coor_xl));
177: h = (int)((XiWin->h)*y*(draw->port_yr - draw->port_yl)/(draw->coor_yr - draw->coor_yl));
178: PetscFree(XiWin->font);
179: PetscDrawXiFontFixed(XiWin,w,h,&XiWin->font);
180: return 0;
181: }
183: static PetscErrorCode PetscDrawStringGetSize_X(PetscDraw draw,PetscReal *x,PetscReal *y)
184: {
185: PetscDraw_X *XiWin = (PetscDraw_X*)draw->data;
186: PetscReal w,h;
188: w = XiWin->font->font_w; h = XiWin->font->font_h;
189: if (x) *x = w*(draw->coor_xr - draw->coor_xl)/((XiWin->w)*(draw->port_xr - draw->port_xl));
190: if (y) *y = h*(draw->coor_yr - draw->coor_yl)/((XiWin->h)*(draw->port_yr - draw->port_yl));
191: return 0;
192: }
194: static PetscErrorCode PetscDrawString_X(PetscDraw draw,PetscReal x,PetscReal y,int c,const char chrs[])
195: {
196: PetscDraw_X *XiWin = (PetscDraw_X*)draw->data;
197: int xx,yy,descent = XiWin->font->font_descent;
198: size_t len;
199: char *substr;
200: PetscToken token;
202: xx = XTRANS(draw,XiWin,x);
203: yy = YTRANS(draw,XiWin,y);
204: PetscDrawXiSetColor(XiWin,c);
206: PetscTokenCreate(chrs,'\n',&token);
207: PetscTokenFind(token,&substr);
208: while (substr) {
209: PetscStrlen(substr,&len);
210: XDrawString(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,xx,yy-descent,substr,len);
211: yy += XiWin->font->font_h;
212: PetscTokenFind(token,&substr);
213: }
214: PetscTokenDestroy(&token);
215: return 0;
216: }
218: static PetscErrorCode PetscDrawStringVertical_X(PetscDraw draw,PetscReal x,PetscReal y,int c,const char text[])
219: {
220: PetscDraw_X *XiWin = (PetscDraw_X*)draw->data;
221: int xx,yy,offset = XiWin->font->font_h - XiWin->font->font_descent;
222: char chr[2] = {0, 0};
224: xx = XTRANS(draw,XiWin,x);
225: yy = YTRANS(draw,XiWin,y);
226: PetscDrawXiSetColor(XiWin,c);
227: while ((chr[0] = *text++)) {
228: XDrawString(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,xx,yy+offset,chr,1);
229: yy += XiWin->font->font_h;
230: }
231: return 0;
232: }
234: static PetscErrorCode PetscDrawFlush_X(PetscDraw draw)
235: {
236: PetscDraw_X *XiWin = (PetscDraw_X*)draw->data;
237: PetscMPIInt rank;
239: /* make sure the X server processed requests from all processes */
240: PetscDrawCollectiveBegin(draw);
241: XSync(XiWin->disp,False);
242: PetscDrawCollectiveEnd(draw);
243: MPI_Barrier(PetscObjectComm((PetscObject)draw));
245: /* transfer pixmap contents to window (only the first process does this) */
246: if (XiWin->drw && XiWin->win) {
247: MPI_Comm_rank(PetscObjectComm((PetscObject)draw),&rank);
248: PetscDrawCollectiveBegin(draw);
249: if (rank == 0) XCopyArea(XiWin->disp,XiWin->drw,XiWin->win,XiWin->gc.set,0,0,XiWin->w,XiWin->h,0,0);
250: if (rank == 0) XSync(XiWin->disp,False);
251: PetscDrawCollectiveEnd(draw);
252: MPI_Barrier(PetscObjectComm((PetscObject)draw));
253: }
254: return 0;
255: }
257: static PetscErrorCode PetscDrawClear_X(PetscDraw draw)
258: {
259: PetscDraw_X *XiWin = (PetscDraw_X*)draw->data;
260: int xmax = XiWin->w-1, ymax = XiWin->h-1;
261: PetscReal xl = draw->port_xl, yl = draw->port_yl;
262: PetscReal xr = draw->port_xr, yr = draw->port_yr;
263: PetscMPIInt rank;
265: /* make sure the X server processed requests from all processes */
266: PetscDrawCollectiveBegin(draw);
267: XSync(XiWin->disp,False);
268: PetscDrawCollectiveEnd(draw);
269: MPI_Barrier(PetscObjectComm((PetscObject)draw));
271: /* only the first process handles the clearing business */
272: PetscDrawCollectiveBegin(draw);
273: MPI_Comm_rank(PetscObjectComm((PetscObject)draw),&rank);
274: if (rank == 0) {
275: int xa = (int)(xl*xmax), ya = ymax - (int)(yr*ymax);
276: int xb = (int)(xr*xmax), yb = ymax - (int)(yl*ymax);
277: unsigned int w = (unsigned int)(xb + 1 - xa);
278: unsigned int h = (unsigned int)(yb + 1 - ya);
279: PetscDrawXiSetPixVal(XiWin,XiWin->background);
280: XFillRectangle(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,xa,ya,w,h);
281: XSync(XiWin->disp,False);
282: }
283: PetscDrawCollectiveEnd(draw);
284: MPI_Barrier(PetscObjectComm((PetscObject)draw));
285: return 0;
286: }
288: static PetscErrorCode PetscDrawSetDoubleBuffer_X(PetscDraw draw)
289: {
290: PetscDraw_X *win = (PetscDraw_X*)draw->data;
291: PetscMPIInt rank;
293: if (win->drw) return 0;
294: MPI_Comm_rank(PetscObjectComm((PetscObject)draw),&rank);
296: PetscDrawCollectiveBegin(draw);
297: if (rank == 0) PetscDrawXiQuickPixmap(win);
298: PetscDrawCollectiveEnd(draw);
299: MPI_Bcast(&win->drw,1,MPI_UNSIGNED_LONG,0,PetscObjectComm((PetscObject)draw));
300: return 0;
301: }
303: static PetscErrorCode PetscDrawGetPopup_X(PetscDraw draw,PetscDraw *popup)
304: {
305: PetscDraw_X *win = (PetscDraw_X*)draw->data;
306: PetscBool flg = PETSC_TRUE;
308: PetscOptionsGetBool(((PetscObject)draw)->options,((PetscObject)draw)->prefix,"-draw_popup",&flg,NULL);
309: if (!flg || !win->win) {*popup = NULL; return 0;}
311: PetscDrawCreate(PetscObjectComm((PetscObject)draw),draw->display,NULL,win->x,win->y+win->h+10,220,220,popup);
312: PetscObjectSetOptionsPrefix((PetscObject)*popup,"popup_");
313: PetscObjectAppendOptionsPrefix((PetscObject)*popup,((PetscObject)draw)->prefix);
314: PetscDrawSetType(*popup,PETSC_DRAW_X);
315: draw->popup = *popup;
316: return 0;
317: }
319: static PetscErrorCode PetscDrawSetTitle_X(PetscDraw draw,const char title[])
320: {
321: PetscDraw_X *win = (PetscDraw_X*)draw->data;
322: PetscMPIInt rank;
324: if (!win->win) return 0;
325: MPI_Comm_rank(PetscObjectComm((PetscObject)draw),&rank);
326: PetscDrawCollectiveBegin(draw);
327: if (rank == 0) {
328: size_t len;
329: XTextProperty prop;
330: PetscStrlen(title,&len);
331: XGetWMName(win->disp,win->win,&prop);
332: XFree((void*)prop.value);
333: prop.value = (unsigned char*)title;
334: prop.nitems = (long)len;
335: XSetWMName(win->disp,win->win,&prop);
336: }
337: PetscDrawCollectiveEnd(draw);
338: return 0;
339: }
341: static PetscErrorCode PetscDrawCheckResizedWindow_X(PetscDraw draw)
342: {
343: PetscDraw_X *win = (PetscDraw_X*)draw->data;
344: int xywh[4];
345: PetscMPIInt rank;
347: if (!win->win) return 0;
348: MPI_Comm_rank(PetscObjectComm((PetscObject)draw),&rank);
350: PetscDrawCollectiveBegin(draw);
351: if (rank == 0) PetscDrawXiGetGeometry(win,xywh,xywh+1,xywh+2,xywh+3);
352: PetscDrawCollectiveEnd(draw);
353: MPI_Bcast(xywh,4,MPI_INT,0,PetscObjectComm((PetscObject)draw));
355: /* record new window position */
356: draw->x = win->x = xywh[0];
357: draw->y = win->y = xywh[1];
358: if (xywh[2] == win->w && xywh[3] == win->h) return 0;
359: /* record new window sizes */
360: draw->w = win->w = xywh[2];
361: draw->h = win->h = xywh[3];
363: /* recreate pixmap (only first processor does this) */
364: PetscDrawCollectiveBegin(draw);
365: if (rank == 0 && win->drw) PetscDrawXiQuickPixmap(win);
366: PetscDrawCollectiveEnd(draw);
367: MPI_Bcast(&win->drw,1,MPI_UNSIGNED_LONG,0,PetscObjectComm((PetscObject)draw));
368: /* reset the clipping */
369: PetscDrawSetViewport_X(draw,draw->port_xl,draw->port_yl,draw->port_xr,draw->port_yr);
370: return 0;
371: }
373: static PetscErrorCode PetscDrawResizeWindow_X(PetscDraw draw,int w,int h)
374: {
375: PetscDraw_X *win = (PetscDraw_X*)draw->data;
376: PetscMPIInt rank;
378: if (w == win->w && h == win->h) return 0;
379: MPI_Comm_rank(PetscObjectComm((PetscObject)draw),&rank);
381: if (win->win) {
382: PetscDrawCollectiveBegin(draw);
383: if (rank == 0) PetscDrawXiResizeWindow(win,w,h);
384: PetscDrawCollectiveEnd(draw);
385: PetscDrawCheckResizedWindow_X(draw);
386: } else if (win->drw) {
387: draw->w = win->w = w; draw->h = win->h = h;
388: /* recreate pixmap (only first processor does this) */
389: PetscDrawCollectiveBegin(draw);
390: if (rank == 0) PetscDrawXiQuickPixmap(win);
391: MPI_Bcast(&win->drw,1,MPI_UNSIGNED_LONG,0,PetscObjectComm((PetscObject)draw));
392: /* reset the clipping */
393: PetscDrawCollectiveEnd(draw);
394: PetscDrawSetViewport_X(draw,draw->port_xl,draw->port_yl,draw->port_xr,draw->port_yr);
395: }
396: return 0;
397: }
399: #include <X11/cursorfont.h>
401: static PetscErrorCode PetscDrawGetMouseButton_X(PetscDraw draw,PetscDrawButton *button,PetscReal *x_user,PetscReal *y_user,PetscReal *x_phys,PetscReal *y_phys)
402: {
403: PetscDraw_X *win = (PetscDraw_X*)draw->data;
404: Cursor cursor;
405: XEvent report;
406: Window root,child;
407: int root_x,root_y,px=0,py=0;
408: unsigned int w,h,border,depth;
409: unsigned int keys_button;
410: PetscMPIInt rank;
411: PetscReal xx,yy;
413: *button = PETSC_BUTTON_NONE;
414: if (!win->win) return 0;
415: MPI_Comm_rank(PetscObjectComm((PetscObject)draw),&rank);
417: PetscDrawCollectiveBegin(draw);
418: if (rank) goto finally;
420: /* change cursor to indicate input */
422: XDefineCursor(win->disp,win->win,cursor);
423: /* wait for mouse button events */
424: XSelectInput(win->disp,win->win,ButtonPressMask|ButtonReleaseMask);
425: while (XCheckTypedEvent(win->disp,ButtonPress,&report));
426: XMaskEvent(win->disp,ButtonReleaseMask,&report);
427: /* get mouse pointer coordinates */
428: XQueryPointer(win->disp,report.xmotion.window,&root,&child,&root_x,&root_y,&px,&py,&keys_button);
429: /* the user may resize the window before pressing the mouse button */
430: XGetGeometry(win->disp,win->win,&root,&root_x,&root_y,&w,&h,&border,&depth);
431: /* cleanup input event handler and cursor */
432: XSelectInput(win->disp,win->win,NoEventMask);
433: XUndefineCursor(win->disp,win->win);
434: XFreeCursor(win->disp, cursor);
435: XSync(win->disp,False);
437: switch (report.xbutton.button) {
438: case Button1: *button = PETSC_BUTTON_LEFT; break;
439: case Button2: *button = PETSC_BUTTON_CENTER; break;
440: case Button3: *button = PETSC_BUTTON_RIGHT; break;
441: case Button4: *button = PETSC_BUTTON_WHEEL_UP; break;
442: case Button5: *button = PETSC_BUTTON_WHEEL_DOWN; break;
443: }
444: if (report.xbutton.state & ShiftMask) {
445: switch (report.xbutton.button) {
446: case Button1: *button = PETSC_BUTTON_LEFT_SHIFT; break;
447: case Button2: *button = PETSC_BUTTON_CENTER_SHIFT; break;
448: case Button3: *button = PETSC_BUTTON_RIGHT_SHIFT; break;
449: }
450: }
451: xx = ((PetscReal)px)/w;
452: yy = 1 - ((PetscReal)py)/h;
453: if (x_user) *x_user = draw->coor_xl + (xx - draw->port_xl)*(draw->coor_xr - draw->coor_xl)/(draw->port_xr - draw->port_xl);
454: if (y_user) *y_user = draw->coor_yl + (yy - draw->port_yl)*(draw->coor_yr - draw->coor_yl)/(draw->port_yr - draw->port_yl);
455: if (x_phys) *x_phys = xx;
456: if (y_phys) *y_phys = yy;
458: finally:
459: PetscDrawCollectiveEnd(draw);
460: PetscDrawCheckResizedWindow_X(draw);
461: return 0;
462: }
464: static PetscErrorCode PetscDrawPause_X(PetscDraw draw)
465: {
466: PetscDraw_X *win = (PetscDraw_X*)draw->data;
468: if (!win->win) return 0;
469: if (draw->pause > 0) PetscSleep(draw->pause);
470: else if (draw->pause == -1) {
471: PetscDrawButton button = PETSC_BUTTON_NONE;
472: PetscDrawGetMouseButton(draw,&button,NULL,NULL,NULL,NULL);
473: if (button == PETSC_BUTTON_CENTER) draw->pause = 0;
474: }
475: return 0;
476: }
478: static PetscErrorCode PetscDrawDestroy_X(PetscDraw draw)
479: {
480: PetscDraw_X *win = (PetscDraw_X*)draw->data;
482: PetscDrawDestroy(&draw->popup);
483: PetscDrawXiClose(win);
484: PetscFree(draw->data);
485: return 0;
486: }
488: static PetscErrorCode PetscDrawGetSingleton_X(PetscDraw,PetscDraw*);
489: static PetscErrorCode PetscDrawRestoreSingleton_X(PetscDraw,PetscDraw*);
490: PETSC_INTERN PetscErrorCode PetscDrawGetImage_X(PetscDraw,unsigned char[][3],unsigned int*,unsigned int*,unsigned char*[]);
492: static struct _PetscDrawOps DvOps = { PetscDrawSetDoubleBuffer_X,
493: PetscDrawFlush_X,
494: PetscDrawLine_X,
495: NULL,
496: NULL,
497: PetscDrawPoint_X,
498: NULL,
499: PetscDrawString_X,
500: PetscDrawStringVertical_X,
501: PetscDrawStringSetSize_X,
502: PetscDrawStringGetSize_X,
503: PetscDrawSetViewport_X,
504: PetscDrawClear_X,
505: PetscDrawRectangle_X,
506: PetscDrawTriangle_X,
507: PetscDrawEllipse_X,
508: PetscDrawGetMouseButton_X,
509: PetscDrawPause_X,
510: NULL,
511: NULL,
512: PetscDrawGetPopup_X,
513: PetscDrawSetTitle_X,
514: PetscDrawCheckResizedWindow_X,
515: PetscDrawResizeWindow_X,
516: PetscDrawDestroy_X,
517: NULL,
518: PetscDrawGetSingleton_X,
519: PetscDrawRestoreSingleton_X,
520: NULL,
521: PetscDrawGetImage_X,
522: NULL,
523: PetscDrawArrow_X,
524: PetscDrawCoordinateToPixel_X,
525: PetscDrawPixelToCoordinate_X,
526: PetscDrawPointPixel_X,
527: NULL};
529: static PetscErrorCode PetscDrawGetSingleton_X(PetscDraw draw,PetscDraw *sdraw)
530: {
531: PetscDraw_X *Xwin = (PetscDraw_X*)draw->data,*sXwin;
533: PetscDrawCreate(PETSC_COMM_SELF,draw->display,draw->title,draw->x,draw->y,draw->w,draw->h,sdraw);
534: PetscObjectChangeTypeName((PetscObject)*sdraw,PETSC_DRAW_X);
535: PetscMemcpy((*sdraw)->ops,&DvOps,sizeof(DvOps));
537: if (draw->popup) PetscDrawGetSingleton(draw->popup,&(*sdraw)->popup);
538: (*sdraw)->pause = draw->pause;
539: (*sdraw)->coor_xl = draw->coor_xl;
540: (*sdraw)->coor_xr = draw->coor_xr;
541: (*sdraw)->coor_yl = draw->coor_yl;
542: (*sdraw)->coor_yr = draw->coor_yr;
543: (*sdraw)->port_xl = draw->port_xl;
544: (*sdraw)->port_xr = draw->port_xr;
545: (*sdraw)->port_yl = draw->port_yl;
546: (*sdraw)->port_yr = draw->port_yr;
548: /* share drawables (windows and/or pixmap) from the parent draw */
549: PetscNewLog(*sdraw,&sXwin);
550: (*sdraw)->data = (void*)sXwin;
551: PetscDrawXiInit(sXwin,draw->display);
552: if (Xwin->win) {
553: PetscDrawXiQuickWindowFromWindow(sXwin,Xwin->win);
554: sXwin->drw = Xwin->drw; /* XXX If the window is ever resized, this is wrong! */
555: } else if (Xwin->drw) {
556: PetscDrawXiColormap(sXwin);
557: sXwin->drw = Xwin->drw;
558: }
559: PetscDrawXiGetGeometry(sXwin,&sXwin->x,&sXwin->y,&sXwin->w,&sXwin->h);
560: (*sdraw)->x = sXwin->x; (*sdraw)->y = sXwin->y;
561: (*sdraw)->w = sXwin->w; (*sdraw)->h = sXwin->h;
562: return 0;
563: }
565: static PetscErrorCode PetscDrawRestoreSingleton_X(PetscDraw draw,PetscDraw *sdraw)
566: {
567: if (draw->popup && (*sdraw)->popup) {
568: PetscBool isdrawx;
569: PetscDraw_X *pXwin = (PetscDraw_X*)draw->popup->data;
570: PetscDraw_X *sXwin = (PetscDraw_X*)(*sdraw)->popup->data;
571: PetscObjectTypeCompare((PetscObject)draw->popup,PETSC_DRAW_X,&isdrawx);
572: if (!isdrawx) goto finally;
573: PetscObjectTypeCompare((PetscObject)(*sdraw)->popup,PETSC_DRAW_X,&isdrawx);
574: if (!isdrawx) goto finally;
575: if (sXwin->win == pXwin->win) {
576: PetscDrawRestoreSingleton(draw->popup,&(*sdraw)->popup);
577: }
578: }
579: finally:
580: PetscDrawDestroy(sdraw);
581: return 0;
582: }
584: static PetscErrorCode PetscDrawXGetDisplaySize_Private(const char name[],int *width,int *height,PetscBool *has_display)
585: {
586: Display *display;
588: display = XOpenDisplay(name);
589: if (!display) {
590: *width = *height = 0;
591: (*PetscErrorPrintf)("Unable to open display on %s\n\
592: Make sure your COMPUTE NODES are authorized to connect\n\
593: to this X server and either your DISPLAY variable\n\
594: is set or you use the -display name option\n",name);
595: *has_display = PETSC_FALSE;
596: return 0;
597: }
598: *has_display = PETSC_TRUE;
599: *width = (int)DisplayWidth(display,DefaultScreen(display));
600: *height = (int)DisplayHeight(display,DefaultScreen(display));
601: XCloseDisplay(display);
602: return 0;
603: }
605: /*MC
606: PETSC_DRAW_X - PETSc graphics device that uses either X windows or its virtual version Xvfb
608: Options Database Keys:
609: + -display <display> - sets the display to use
610: . -x_virtual - forces use of a X virtual display Xvfb that will not display anything but -draw_save will still work.
611: Xvfb is automatically started up in PetscSetDisplay() with this option
612: . -draw_size w,h - percentage of screen (either 1, .5, .3, .25), or size in pixels
613: . -geometry x,y,w,h - set location and size in pixels
614: . -draw_virtual - do not open a window (draw on a pixmap), -draw_save will still work
615: - -draw_double_buffer - avoid window flickering (draw on pixmap and flush to window)
617: Level: beginner
619: .seealso: `PetscDrawOpenX()`, `PetscDrawSetDisplay()`, `PetscDrawSetFromOptions()`
621: M*/
623: PETSC_EXTERN PetscErrorCode PetscDrawCreate_X(PetscDraw draw)
624: {
625: PetscDraw_X *Xwin;
626: PetscMPIInt rank;
627: int x = draw->x,y = draw->y,w = draw->w,h = draw->h;
628: static int xavailable = 0,yavailable = 0,ybottom = 0,xmax = 0,ymax = 0;
629: PetscBool set,dvirtual = PETSC_FALSE,doublebuffer = PETSC_TRUE,has_display;
630: PetscInt xywh[4],osize = 4,nsizes=2;
631: PetscReal sizes[2] = {.3,.3};
632: static size_t DISPLAY_LENGTH = 265;
634: /* get the display variable */
635: if (!draw->display) {
636: PetscMalloc1(DISPLAY_LENGTH,&draw->display);
637: PetscGetDisplay(draw->display,DISPLAY_LENGTH);
638: }
640: /* initialize the display size */
641: if (!xmax) {
642: PetscDrawXGetDisplaySize_Private(draw->display,&xmax,&ymax,&has_display);
643: /* if some processors fail on this and others succed then this is a problem ! */
644: if (!has_display) {
645: (*PetscErrorPrintf)("PETSc unable to use X windows\nproceeding without graphics\n");
646: PetscDrawSetType(draw,PETSC_DRAW_NULL);
647: return 0;
648: }
649: }
651: /* allow user to set size of drawable */
652: PetscOptionsGetRealArray(((PetscObject)draw)->options,((PetscObject)draw)->prefix,"-draw_size",sizes,&nsizes,&set);
653: if (set && nsizes == 1 && sizes[0] > 1.0) sizes[1] = sizes[0];
654: if (set) {
655: if (sizes[0] > 1.0) w = (int)sizes[0];
656: else if (sizes[0] == 1.0) w = PETSC_DRAW_FULL_SIZE;
657: else if (sizes[0] == .5) w = PETSC_DRAW_HALF_SIZE;
658: else if (sizes[0] == .3) w = PETSC_DRAW_THIRD_SIZE;
659: else if (sizes[0] == .25) w = PETSC_DRAW_QUARTER_SIZE;
660: if (sizes[1] > 1.0) h = (int)sizes[1];
661: else if (sizes[1] == 1.0) h = PETSC_DRAW_FULL_SIZE;
662: else if (sizes[1] == .5) h = PETSC_DRAW_HALF_SIZE;
663: else if (sizes[1] == .3) h = PETSC_DRAW_THIRD_SIZE;
664: else if (sizes[1] == .25) h = PETSC_DRAW_QUARTER_SIZE;
665: }
666: if (w == PETSC_DECIDE || w == PETSC_DEFAULT) w = draw->w = 300;
667: if (h == PETSC_DECIDE || h == PETSC_DEFAULT) h = draw->h = 300;
668: switch (w) {
669: case PETSC_DRAW_FULL_SIZE: w = draw->w = (xmax - 10); break;
670: case PETSC_DRAW_HALF_SIZE: w = draw->w = (xmax - 20)/2; break;
671: case PETSC_DRAW_THIRD_SIZE: w = draw->w = (xmax - 30)/3; break;
672: case PETSC_DRAW_QUARTER_SIZE: w = draw->w = (xmax - 40)/4; break;
673: }
674: switch (h) {
675: case PETSC_DRAW_FULL_SIZE: h = draw->h = (ymax - 10); break;
676: case PETSC_DRAW_HALF_SIZE: h = draw->h = (ymax - 20)/2; break;
677: case PETSC_DRAW_THIRD_SIZE: h = draw->h = (ymax - 30)/3; break;
678: case PETSC_DRAW_QUARTER_SIZE: h = draw->h = (ymax - 40)/4; break;
679: }
681: PetscOptionsGetBool(((PetscObject)draw)->options,((PetscObject)draw)->prefix,"-draw_virtual",&dvirtual,NULL);
683: if (!dvirtual) {
685: /* allow user to set location and size of window */
686: xywh[0] = x; xywh[1] = y; xywh[2] = w; xywh[3] = h;
687: PetscOptionsGetIntArray(((PetscObject)draw)->options,((PetscObject)draw)->prefix,"-geometry",xywh,&osize,NULL);
688: x = (int)xywh[0]; y = (int)xywh[1]; w = (int)xywh[2]; h = (int)xywh[3];
689: if (w == PETSC_DECIDE || w == PETSC_DEFAULT) w = 300;
690: if (h == PETSC_DECIDE || h == PETSC_DEFAULT) h = 300;
691: draw->x = x; draw->y = y; draw->w = w; draw->h = h;
693: if (draw->x == PETSC_DECIDE || draw->y == PETSC_DECIDE) {
694: /*
695: PETSc tries to place windows starting in the upper left corner
696: and moving across to the right.
698: +0,0-------------------------------------------+
699: | Region used so far +xavailable,yavailable |
700: | | |
701: | | |
702: +--------------------- +ybottom |
703: | |
704: | |
705: +----------------------------------------------+xmax,ymax
707: */
708: /* First: can we add it to the right? */
709: if (xavailable + w + 10 <= xmax) {
710: x = xavailable;
711: y = yavailable;
712: ybottom = PetscMax(ybottom,y + h + 30);
713: } else {
714: /* No, so add it below on the left */
715: xavailable = x = 0;
716: yavailable = y = ybottom;
717: ybottom = ybottom + h + 30;
718: }
719: }
720: /* update available region */
721: xavailable = PetscMax(xavailable,x + w + 10);
722: if (xavailable >= xmax) {
723: xavailable = 0;
724: yavailable = yavailable + h + 30;
725: ybottom = yavailable;
726: }
727: if (yavailable >= ymax) {
728: y = 0;
729: yavailable = 0;
730: ybottom = 0;
731: }
733: } /* endif (!dvirtual) */
735: MPI_Comm_rank(PetscObjectComm((PetscObject)draw),&rank);
738: PetscNewLog(draw,&Xwin);
739: PetscMemcpy(draw->ops,&DvOps,sizeof(DvOps));
740: draw->data = (void*)Xwin;
742: PetscDrawXiInit(Xwin,draw->display);
743: if (!dvirtual) {
744: Xwin->x = x; Xwin->y = y;
745: Xwin->w = w; Xwin->h = h;
746: if (rank == 0) PetscDrawXiQuickWindow(Xwin,draw->title,x,y,w,h);
747: MPI_Bcast(&Xwin->win,1,MPI_UNSIGNED_LONG,0,PetscObjectComm((PetscObject)draw));
748: if (rank) PetscDrawXiQuickWindowFromWindow(Xwin,Xwin->win);
749: } else {
750: Xwin->x = 0; Xwin->y = 0;
751: Xwin->w = w; Xwin->h = h;
752: PetscDrawXiColormap(Xwin);
753: if (rank == 0) PetscDrawXiQuickPixmap(Xwin);
754: MPI_Bcast(&Xwin->drw,1,MPI_UNSIGNED_LONG,0,PetscObjectComm((PetscObject)draw));
755: }
756: PetscDrawXiGetGeometry(Xwin,&Xwin->x,&Xwin->y,&Xwin->w,&Xwin->h);
757: draw->x = Xwin->x; draw->y = Xwin->y;
758: draw->w = Xwin->w; draw->h = Xwin->h;
760: PetscOptionsGetBool(((PetscObject)draw)->options,((PetscObject)draw)->prefix,"-draw_double_buffer",&doublebuffer,NULL);
761: if (doublebuffer) PetscDrawSetDoubleBuffer(draw);
762: return 0;
763: }