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