Actual source code: xops.c
2: /*
3: Defines the operations for the X PetscDraw implementation.
4: */
6: #include ../src/sys/draw/impls/x/ximpl.h
8: /*
9: These macros transform from the users coordinates to the
10: X-window pixel coordinates.
11: */
12: #define XTRANS(draw,xwin,x) \
13: (int)(((xwin)->w)*((draw)->port_xl + (((x - (draw)->coor_xl)*\
14: ((draw)->port_xr - (draw)->port_xl))/\
15: ((draw)->coor_xr - (draw)->coor_xl))))
16: #define YTRANS(draw,xwin,y) \
17: (int)(((xwin)->h)*(1.0-(draw)->port_yl - (((y - (draw)->coor_yl)*\
18: ((draw)->port_yr - (draw)->port_yl))/\
19: ((draw)->coor_yr - (draw)->coor_yl))))
23: PetscErrorCode PetscDrawLine_X(PetscDraw draw,PetscReal xl,PetscReal yl,PetscReal xr,PetscReal yr,int cl)
24: {
25: PetscDraw_X* XiWin = (PetscDraw_X*)draw->data;
26: int x1,y_1,x2,y2;
29: XiSetColor(XiWin,cl);
30: x1 = XTRANS(draw,XiWin,xl); x2 = XTRANS(draw,XiWin,xr);
31: y_1 = YTRANS(draw,XiWin,yl); y2 = YTRANS(draw,XiWin,yr);
32: XDrawLine(XiWin->disp,XiDrawable(XiWin),XiWin->gc.set,x1,y_1,x2,y2);
33: return(0);
34: }
38: static PetscErrorCode PetscDrawPoint_X(PetscDraw draw,PetscReal x,PetscReal y,int c)
39: {
40: int xx,yy;
41: PetscDraw_X* XiWin = (PetscDraw_X*)draw->data;
44: xx = XTRANS(draw,XiWin,x); yy = YTRANS(draw,XiWin,y);
45: XiSetColor(XiWin,c);
46: XDrawPoint(XiWin->disp,XiDrawable(XiWin),XiWin->gc.set,xx,yy);
47: return(0);
48: }
52: static PetscErrorCode PetscDrawRectangle_X(PetscDraw draw,PetscReal xl,PetscReal yl,PetscReal xr,PetscReal yr,int c1,int c2,int c3,int c4)
53: {
54: PetscDraw_X* XiWin = (PetscDraw_X*)draw->data;
55: int x1,y_1,w,h,c = (c1 + c2 + c3 + c4)/4;
58: XiSetColor(XiWin,c);
59: x1 = XTRANS(draw,XiWin,xl); w = XTRANS(draw,XiWin,xr) - x1;
60: y_1 = YTRANS(draw,XiWin,yr); h = YTRANS(draw,XiWin,yl) - y_1;
61: if (w <= 0) w = 1; if (h <= 0) h = 1;
62: XFillRectangle(XiWin->disp,XiDrawable(XiWin),XiWin->gc.set,x1,y_1,w,h);
63: return(0);
64: }
68: static PetscErrorCode PetscDrawEllipse_X(PetscDraw Win, PetscReal x, PetscReal y, PetscReal a, PetscReal b, int c)
69: {
70: PetscDraw_X* XiWin = (PetscDraw_X*) Win->data;
71: int xA, yA, w, h;
74: XiSetColor(XiWin, c);
75: xA = XTRANS(Win, XiWin, x - a/2.0); w = XTRANS(Win, XiWin, x + a/2.0) - xA;
76: yA = YTRANS(Win, XiWin, y + b/2.0); h = YTRANS(Win, XiWin, y - b/2.0) - yA;
77: XFillArc(XiWin->disp, XiDrawable(XiWin), XiWin->gc.set, xA, yA, w, h, 0, 23040);
78: return(0);
79: }
81: EXTERN PetscErrorCode PetscDrawInterpolatedTriangle_X(PetscDraw_X*,int,int,int,int,int,int,int,int,int);
85: static PetscErrorCode PetscDrawTriangle_X(PetscDraw draw,PetscReal X1,PetscReal Y_1,PetscReal X2,
86: PetscReal Y2,PetscReal X3,PetscReal Y3,int c1,int c2,int c3)
87: {
88: PetscDraw_X* XiWin = (PetscDraw_X*)draw->data;
92: if (c1 == c2 && c2 == c3) {
93: XPoint pt[3];
94: XiSetColor(XiWin,c1);
95: pt[0].x = XTRANS(draw,XiWin,X1);
96: pt[0].y = YTRANS(draw,XiWin,Y_1);
97: pt[1].x = XTRANS(draw,XiWin,X2);
98: pt[1].y = YTRANS(draw,XiWin,Y2);
99: pt[2].x = XTRANS(draw,XiWin,X3);
100: pt[2].y = YTRANS(draw,XiWin,Y3);
101: XFillPolygon(XiWin->disp,XiDrawable(XiWin),XiWin->gc.set,pt,3,Convex,CoordModeOrigin);
102: } else {
103: int x1,y_1,x2,y2,x3,y3;
104: x1 = XTRANS(draw,XiWin,X1);
105: y_1 = YTRANS(draw,XiWin,Y_1);
106: x2 = XTRANS(draw,XiWin,X2);
107: y2 = YTRANS(draw,XiWin,Y2);
108: x3 = XTRANS(draw,XiWin,X3);
109: y3 = YTRANS(draw,XiWin,Y3);
110: PetscDrawInterpolatedTriangle_X(XiWin,x1,y_1,c1,x2,y2,c2,x3,y3,c3);
111: }
112: return(0);
113: }
117: static PetscErrorCode PetscDrawString_X(PetscDraw draw,PetscReal x,PetscReal y,int c,const char chrs[])
118: {
120: int xx,yy;
121: size_t len;
122: PetscDraw_X *XiWin = (PetscDraw_X*)draw->data;
123: char *substr;
124: PetscToken token;
127: xx = XTRANS(draw,XiWin,x); yy = YTRANS(draw,XiWin,y);
128: XiSetColor(XiWin,c);
129:
130: PetscTokenCreate(chrs,'\n',&token);
131: PetscTokenFind(token,&substr);
132: PetscStrlen(substr,&len);
133: XDrawString(XiWin->disp,XiDrawable(XiWin),XiWin->gc.set,xx,yy - XiWin->font->font_descent,substr,len);
134: PetscTokenFind(token,&substr);
135: while (substr) {
136: yy += 4*XiWin->font->font_descent;
137: PetscStrlen(substr,&len);
138: XDrawString(XiWin->disp,XiDrawable(XiWin),XiWin->gc.set,xx,yy - XiWin->font->font_descent,substr,len);
139: PetscTokenFind(token,&substr);
140: }
141: PetscTokenDestroy(token);
143: return(0);
144: }
146: EXTERN PetscErrorCode XiFontFixed(PetscDraw_X*,int,int,XiFont **);
150: static PetscErrorCode PetscDrawStringSetSize_X(PetscDraw draw,PetscReal x,PetscReal y)
151: {
152: PetscDraw_X* XiWin = (PetscDraw_X*)draw->data;
154: int w,h;
157: w = (int)((XiWin->w)*x*(draw->port_xr - draw->port_xl)/(draw->coor_xr - draw->coor_xl));
158: h = (int)((XiWin->h)*y*(draw->port_yr - draw->port_yl)/(draw->coor_yr - draw->coor_yl));
159: PetscFree(XiWin->font);
160: XiFontFixed(XiWin,w,h,&XiWin->font);
161: return(0);
162: }
166: PetscErrorCode PetscDrawStringGetSize_X(PetscDraw draw,PetscReal *x,PetscReal *y)
167: {
168: PetscDraw_X* XiWin = (PetscDraw_X*)draw->data;
169: PetscReal w,h;
172: w = XiWin->font->font_w; h = XiWin->font->font_h;
173: *x = w*(draw->coor_xr - draw->coor_xl)/(XiWin->w)*(draw->port_xr - draw->port_xl);
174: *y = h*(draw->coor_yr - draw->coor_yl)/(XiWin->h)*(draw->port_yr - draw->port_yl);
175: return(0);
176: }
180: PetscErrorCode PetscDrawStringVertical_X(PetscDraw draw,PetscReal x,PetscReal y,int c,const char chrs[])
181: {
183: int xx,yy;
184: PetscDraw_X *XiWin = (PetscDraw_X*)draw->data;
185: char tmp[2];
186: PetscReal tw,th;
187: size_t i,n;
188:
190: PetscStrlen(chrs,&n);
191: tmp[1] = 0;
192: XiSetColor(XiWin,c);
193: PetscDrawStringGetSize_X(draw,&tw,&th);
194: xx = XTRANS(draw,XiWin,x);
195: for (i=0; i<n; i++) {
196: tmp[0] = chrs[i];
197: yy = YTRANS(draw,XiWin,y-th*i);
198: XDrawString(XiWin->disp,XiDrawable(XiWin),XiWin->gc.set,
199: xx,yy - XiWin->font->font_descent,tmp,1);
200: }
201: return(0);
202: }
206: static PetscErrorCode PetscDrawFlush_X(PetscDraw draw)
207: {
208: PetscDraw_X* XiWin = (PetscDraw_X*)draw->data;
211: if (XiWin->drw) {
212: XCopyArea(XiWin->disp,XiWin->drw,XiWin->win,XiWin->gc.set,0,0,XiWin->w,XiWin->h,0,0);
213: }
214: XFlush(XiWin->disp); XSync(XiWin->disp,False);
215: return(0);
216: }
220: static PetscErrorCode PetscDrawSynchronizedFlush_X(PetscDraw draw)
221: {
223: PetscMPIInt rank;
224: PetscDraw_X* XiWin = (PetscDraw_X*)draw->data;
227: XFlush(XiWin->disp);
228: if (XiWin->drw) {
229: MPI_Comm_rank(((PetscObject)draw)->comm,&rank);
230: /* make sure data has actually arrived at server */
231: XSync(XiWin->disp,False);
232: MPI_Barrier(((PetscObject)draw)->comm);
233: if (!rank) {
234: XCopyArea(XiWin->disp,XiWin->drw,XiWin->win,XiWin->gc.set,0,0,XiWin->w,XiWin->h,0,0);
235: XFlush(XiWin->disp);
236: }
237: XSync(XiWin->disp,False);
238: MPI_Barrier(((PetscObject)draw)->comm);
239: } else {
240: MPI_Barrier(((PetscObject)draw)->comm);
241: XSync(XiWin->disp,False);
242: MPI_Barrier(((PetscObject)draw)->comm);
243: }
244: return(0);
245: }
249: static PetscErrorCode PetscDrawSetViewport_X(PetscDraw draw,PetscReal xl,PetscReal yl,PetscReal xr,PetscReal yr)
250: {
251: PetscDraw_X* XiWin = (PetscDraw_X*)draw->data;
252: XRectangle box;
255: box.x = (int)(xl*XiWin->w); box.y = (int)((1.0-yr)*XiWin->h);
256: box.width = (int)((xr-xl)*XiWin->w);box.height = (int)((yr-yl)*XiWin->h);
257: XSetClipRectangles(XiWin->disp,XiWin->gc.set,0,0,&box,1,Unsorted);
258: return(0);
259: }
263: static PetscErrorCode PetscDrawClear_X(PetscDraw draw)
264: {
265: PetscDraw_X* XiWin = (PetscDraw_X*)draw->data;
266: int x, y, w, h;
269: x = (int)(draw->port_xl*XiWin->w);
270: w = (int)((draw->port_xr - draw->port_xl)*XiWin->w);
271: y = (int)((1.0-draw->port_yr)*XiWin->h);
272: h = (int)((draw->port_yr - draw->port_yl)*XiWin->h);
273: XiSetPixVal(XiWin,XiWin->background);
274: XFillRectangle(XiWin->disp,XiDrawable(XiWin),XiWin->gc.set,x,y,w,h);
275: return(0);
276: }
280: static PetscErrorCode PetscDrawSynchronizedClear_X(PetscDraw draw)
281: {
283: PetscMPIInt rank;
284: PetscDraw_X* XiWin = (PetscDraw_X*)draw->data;
287: MPI_Barrier(((PetscObject)draw)->comm);
288: MPI_Comm_rank(((PetscObject)draw)->comm,&rank);
289: if (!rank) {
290: PetscDrawClear_X(draw);
291: }
292: XFlush(XiWin->disp);
293: MPI_Barrier(((PetscObject)draw)->comm);
294: XSync(XiWin->disp,False);
295: MPI_Barrier(((PetscObject)draw)->comm);
296: return(0);
297: }
301: static PetscErrorCode PetscDrawSetDoubleBuffer_X(PetscDraw draw)
302: {
303: PetscDraw_X* win = (PetscDraw_X*)draw->data;
305: PetscMPIInt rank;
308: if (win->drw) return(0);
310: MPI_Comm_rank(((PetscObject)draw)->comm,&rank);
311: if (!rank) {
312: win->drw = XCreatePixmap(win->disp,win->win,win->w,win->h,win->depth);
313: }
314: /* try to make sure it is actually done before passing info to all */
315: XSync(win->disp,False);
316: MPI_Bcast(&win->drw,1,MPI_UNSIGNED_LONG,0,((PetscObject)draw)->comm);
317: return(0);
318: }
320: #include <X11/cursorfont.h>
324: static PetscErrorCode PetscDrawGetMouseButton_X(PetscDraw draw,PetscDrawButton *button,PetscReal* x_user,PetscReal *y_user,PetscReal *x_phys,PetscReal *y_phys)
325: {
326: XEvent report;
327: PetscDraw_X* win = (PetscDraw_X*)draw->data;
328: Window root,child;
329: int root_x,root_y,px,py;
330: unsigned int keys_button;
331: Cursor cursor = 0;
334: /* change cursor to indicate input */
335: if (!cursor) {
336: cursor = XCreateFontCursor(win->disp,XC_hand2);
337: if (!cursor) SETERRQ(PETSC_ERR_LIB,"Unable to create X cursor");
338: }
339: XDefineCursor(win->disp,win->win,cursor);
340: XSelectInput(win->disp,win->win,ButtonPressMask | ButtonReleaseMask);
342: while (XCheckTypedEvent(win->disp,ButtonPress,&report));
343: XMaskEvent(win->disp,ButtonReleaseMask,&report);
344: switch (report.xbutton.button) {
345: case Button1:
346: if (report.xbutton.state & ShiftMask)
347: *button = BUTTON_LEFT_SHIFT;
348: else
349: *button = BUTTON_LEFT;
350: break;
351: case Button2:
352: if (report.xbutton.state & ShiftMask)
353: *button = BUTTON_CENTER_SHIFT;
354: else
355: *button = BUTTON_CENTER;
356: break;
357: case Button3:
358: if (report.xbutton.state & ShiftMask)
359: *button = BUTTON_RIGHT_SHIFT;
360: else
361: *button = BUTTON_RIGHT;
362: break;
363: }
364: XQueryPointer(win->disp,report.xmotion.window,&root,&child,&root_x,&root_y,&px,&py,&keys_button);
366: if (x_phys) *x_phys = ((double)px)/((double)win->w);
367: if (y_phys) *y_phys = 1.0 - ((double)py)/((double)win->h);
369: if (x_user) *x_user = draw->coor_xl + ((((double)px)/((double)win->w)-draw->port_xl))*(draw->coor_xr - draw->coor_xl)/(draw->port_xr - draw->port_xl);
370: if (y_user) *y_user = draw->coor_yl + ((1.0 - ((double)py)/((double)win->h)-draw->port_yl))*(draw->coor_yr - draw->coor_yl)/(draw->port_yr - draw->port_yl);
372: XUndefineCursor(win->disp,win->win);
373: XFlush(win->disp); XSync(win->disp,False);
374: return(0);
375: }
379: static PetscErrorCode PetscDrawPause_X(PetscDraw draw)
380: {
384: if (draw->pause > 0) PetscSleep(draw->pause);
385: else if (draw->pause < 0) {
386: PetscDrawButton button;
387: PetscMPIInt rank;
388: MPI_Comm_rank(((PetscObject)draw)->comm,&rank);
389: if (!rank) {
390: PetscDrawGetMouseButton(draw,&button,0,0,0,0);
391: if (button == BUTTON_CENTER) draw->pause = 0;
392: }
393: MPI_Bcast(&draw->pause,1,MPI_INT,0,((PetscObject)draw)->comm);
394: }
395: return(0);
396: }
400: static PetscErrorCode PetscDrawGetPopup_X(PetscDraw draw,PetscDraw *popup)
401: {
403: PetscDraw_X* win = (PetscDraw_X*)draw->data;
406: PetscDrawOpenX(((PetscObject)draw)->comm,PETSC_NULL,PETSC_NULL,win->x,win->y+win->h+36,150,220,popup);
407: draw->popup = *popup;
408: return(0);
409: }
413: static PetscErrorCode PetscDrawSetTitle_X(PetscDraw draw,const char title[])
414: {
415: PetscDraw_X *win = (PetscDraw_X*)draw->data;
416: XTextProperty prop;
418: size_t len;
421: XGetWMName(win->disp,win->win,&prop);
422: XFree((void*)prop.value);
423: prop.value = (unsigned char *)title;
424: PetscStrlen(title,&len);
425: prop.nitems = (long) len;
426: XSetWMName(win->disp,win->win,&prop);
427: return(0);
428: }
432: static PetscErrorCode PetscDrawResizeWindow_X(PetscDraw draw,int w,int h)
433: {
434: PetscDraw_X *win = (PetscDraw_X*)draw->data;
435: unsigned int ww,hh,border,depth;
436: int x,y;
438: Window root;
441: XResizeWindow(win->disp,win->win,w,h);
442: XGetGeometry(win->disp,win->win,&root,&x,&y,&ww,&hh,&border,&depth);
443: PetscDrawCheckResizedWindow(draw);
444: return(0);
445: }
449: static PetscErrorCode PetscDrawCheckResizedWindow_X(PetscDraw draw)
450: {
451: PetscDraw_X *win = (PetscDraw_X*)draw->data;
453: int x,y;
454: PetscMPIInt rank;
455: Window root;
456: unsigned int w,h,border,depth,geo[2];
457: PetscReal xl,xr,yl,yr;
458: XRectangle box;
461: MPI_Comm_rank(((PetscObject)draw)->comm,&rank);
462: if (!rank) {
463: XSync(win->disp,False);
464: XGetGeometry(win->disp,win->win,&root,&x,&y,geo,geo+1,&border,&depth);
465: }
466: MPI_Bcast(geo,2,MPI_INT,0,((PetscObject)draw)->comm);
467: w = geo[0];
468: h = geo[1];
469: if (w == (unsigned int) win->w && h == (unsigned int) win->h) return(0);
471: /* record new window sizes */
473: win->h = h; win->w = w;
475: /* Free buffer space and create new version (only first processor does this) */
476: if (win->drw) {
477: win->drw = XCreatePixmap(win->disp,win->win,win->w,win->h,win->depth);
478: }
479: /* reset the clipping */
480: xl = draw->port_xl; yl = draw->port_yl;
481: xr = draw->port_xr; yr = draw->port_yr;
482: box.x = (int)(xl*win->w); box.y = (int)((1.0-yr)*win->h);
483: box.width = (int)((xr-xl)*win->w);box.height = (int)((yr-yl)*win->h);
484: XSetClipRectangles(win->disp,win->gc.set,0,0,&box,1,Unsorted);
486: /* try to make sure it is actually done before passing info to all */
487: XSync(win->disp,False);
488: MPI_Bcast(&win->drw,1,MPI_UNSIGNED_LONG,0,((PetscObject)draw)->comm);
489: return(0);
490: }
492: static PetscErrorCode PetscDrawGetSingleton_X(PetscDraw,PetscDraw*);
493: static PetscErrorCode PetscDrawRestoreSingleton_X(PetscDraw,PetscDraw*);
497: PetscErrorCode PetscDrawDestroy_X(PetscDraw draw)
498: {
499: PetscDraw_X *win = (PetscDraw_X*)draw->data;
503: XFreeGC(win->disp,win->gc.set);
504: XCloseDisplay(win->disp);
505: if (draw->popup) {PetscDrawDestroy(draw->popup);}
506: PetscFree(win->font);
507: PetscFree(win);
508: return(0);
509: }
511: static struct _PetscDrawOps DvOps = { PetscDrawSetDoubleBuffer_X,
512: PetscDrawFlush_X,PetscDrawLine_X,
513: 0,
514: 0,
515: PetscDrawPoint_X,
516: 0,
517: PetscDrawString_X,
518: PetscDrawStringVertical_X,
519: PetscDrawStringSetSize_X,
520: PetscDrawStringGetSize_X,
521: PetscDrawSetViewport_X,
522: PetscDrawClear_X,
523: PetscDrawSynchronizedFlush_X,
524: PetscDrawRectangle_X,
525: PetscDrawTriangle_X,
526: PetscDrawEllipse_X,
527: PetscDrawGetMouseButton_X,
528: PetscDrawPause_X,
529: PetscDrawSynchronizedClear_X,
530: 0,
531: 0,
532: PetscDrawGetPopup_X,
533: PetscDrawSetTitle_X,
534: PetscDrawCheckResizedWindow_X,
535: PetscDrawResizeWindow_X,
536: PetscDrawDestroy_X,
537: 0,
538: PetscDrawGetSingleton_X,
539: PetscDrawRestoreSingleton_X };
542: EXTERN PetscErrorCode XiQuickWindow(PetscDraw_X*,char*,char*,int,int,int,int);
543: EXTERN PetscErrorCode XiQuickWindowFromWindow(PetscDraw_X*,char*,Window);
547: static PetscErrorCode PetscDrawGetSingleton_X(PetscDraw draw,PetscDraw *sdraw)
548: {
550: PetscDraw_X *Xwin = (PetscDraw_X*)draw->data,*sXwin;
554: PetscDrawCreate(PETSC_COMM_SELF,draw->display,draw->title,draw->x,draw->y,draw->w,draw->h,sdraw);
555: PetscObjectChangeTypeName((PetscObject)*sdraw,PETSC_DRAW_X);
556: PetscMemcpy((*sdraw)->ops,&DvOps,sizeof(DvOps));
557: (*sdraw)->ops->destroy = 0;
559: (*sdraw)->pause = draw->pause;
560: (*sdraw)->coor_xl = draw->coor_xl;
561: (*sdraw)->coor_xr = draw->coor_xr;
562: (*sdraw)->coor_yl = draw->coor_yl;
563: (*sdraw)->coor_yr = draw->coor_yr;
564: (*sdraw)->port_xl = draw->port_xl;
565: (*sdraw)->port_xr = draw->port_xr;
566: (*sdraw)->port_yl = draw->port_yl;
567: (*sdraw)->port_yr = draw->port_yr;
568: (*sdraw)->popup = draw->popup;
570: /* actually create and open the window */
571: PetscNew(PetscDraw_X,&sXwin);
572: XiQuickWindowFromWindow(sXwin,draw->display,Xwin->win);
574: sXwin->x = Xwin->x;
575: sXwin->y = Xwin->y;
576: sXwin->w = Xwin->w;
577: sXwin->h = Xwin->h;
578: (*sdraw)->data = (void*)sXwin;
579: return(0);
580: }
584: static PetscErrorCode PetscDrawRestoreSingleton_X(PetscDraw draw,PetscDraw *sdraw)
585: {
587: PetscDraw_X *sXwin = (PetscDraw_X*)(*sdraw)->data;
590: XFreeGC(sXwin->disp,sXwin->gc.set);
591: XCloseDisplay(sXwin->disp);
592: if ((*sdraw)->popup) {PetscDrawDestroy((*sdraw)->popup);}
593: PetscStrfree((*sdraw)->title);
594: PetscStrfree((*sdraw)->display);
595: PetscFree(sXwin->font);
596: PetscFree(sXwin);
597: PetscHeaderDestroy(*sdraw);
598: return(0);
599: }
603: PetscErrorCode PetscDrawXGetDisplaySize_Private(const char name[],int *width,int *height)
604: {
605: Display *display;
608: display = XOpenDisplay(name);
609: if (!display) {
610: *width = 0;
611: *height = 0;
612: SETERRQ1(PETSC_ERR_LIB,"Unable to open display on %s\n. Make sure your COMPUTE NODES are authorized to connect \n\
613: to this X server and either your DISPLAY variable\n\
614: is set or you use the -display name option\n",name);
615: }
616: *width = DisplayWidth(display,0);
617: *height = DisplayHeight(display,0);
618: XCloseDisplay(display);
619: return(0);
620: }
625: PetscErrorCode PetscDrawCreate_X(PetscDraw draw)
626: {
627: PetscDraw_X *Xwin;
629: PetscMPIInt rank;
630: PetscInt xywh[4],osize = 4;
631: int x = draw->x,y = draw->y,w = draw->w,h = draw->h;
632: static int xavailable = 0,yavailable = 0,xmax = 0,ymax = 0,ybottom = 0;
633: PetscTruth flg = PETSC_FALSE;
636: if (!draw->display) {
637: PetscMalloc(128*sizeof(char),&draw->display);
638: PetscGetDisplay(draw->display,128);
639: }
641: /*
642: Initialize the display size
643: */
644: if (!xmax) {
645: PetscDrawXGetDisplaySize_Private(draw->display,&xmax,&ymax);
646: /* if some processors fail on this and others succed then this is a problem ! */
647: if (ierr) {
648: (*PetscErrorPrintf)("PETSc unable to use X windows\nproceeding without graphics\n");
649: PetscDrawSetType(draw,PETSC_DRAW_NULL);
650: return(0);
651: }
652: }
654: if (w == PETSC_DECIDE) w = draw->w = 300;
655: if (h == PETSC_DECIDE) h = draw->h = 300;
656: switch (w) {
657: case PETSC_DRAW_FULL_SIZE: w = draw->w = xmax - 10;
658: break;
659: case PETSC_DRAW_HALF_SIZE: w = draw->w = (xmax - 20)/2;
660: break;
661: case PETSC_DRAW_THIRD_SIZE: w = draw->w = (xmax - 30)/3;
662: break;
663: case PETSC_DRAW_QUARTER_SIZE: w = draw->w = (xmax - 40)/4;
664: break;
665: }
666: switch (h) {
667: case PETSC_DRAW_FULL_SIZE: h = draw->h = ymax - 10;
668: break;
669: case PETSC_DRAW_HALF_SIZE: h = draw->h = (ymax - 20)/2;
670: break;
671: case PETSC_DRAW_THIRD_SIZE: h = draw->h = (ymax - 30)/3;
672: break;
673: case PETSC_DRAW_QUARTER_SIZE: h = draw->h = (ymax - 40)/4;
674: break;
675: }
677: /* allow user to set location and size of window */
678: xywh[0] = x; xywh[1] = y; xywh[2] = w; xywh[3] = h;
679: PetscOptionsGetIntArray(PETSC_NULL,"-geometry",xywh,&osize,PETSC_NULL);
680: x = (int) xywh[0]; y = (int) xywh[1]; w = (int) xywh[2]; h = (int) xywh[3];
683: if (draw->x == PETSC_DECIDE || draw->y == PETSC_DECIDE) {
684: /*
685: PETSc tries to place windows starting in the upper left corner and
686: moving across to the right.
687:
688: --------------------------------------------
689: | Region used so far +xavailable,yavailable |
690: | + |
691: | + |
692: |++++++++++++++++++++++ybottom |
693: | |
694: | |
695: |--------------------------------------------|
696: */
697: /* First: can we add it to the right? */
698: if (xavailable+w+10 <= xmax) {
699: x = xavailable;
700: y = yavailable;
701: ybottom = PetscMax(ybottom,y + h + 30);
702: } else {
703: /* No, so add it below on the left */
704: xavailable = 0;
705: x = 0;
706: yavailable = ybottom;
707: y = ybottom;
708: ybottom = ybottom + h + 30;
709: }
710: }
711: /* update available region */
712: xavailable = PetscMax(xavailable,x + w + 10);
713: if (xavailable >= xmax) {
714: xavailable = 0;
715: yavailable = yavailable + h + 30;
716: ybottom = yavailable;
717: }
718: if (yavailable >= ymax) {
719: y = 0;
720: yavailable = 0;
721: ybottom = 0;
722: }
724: PetscMemcpy(draw->ops,&DvOps,sizeof(DvOps));
726: /* actually create and open the window */
727: PetscNew(PetscDraw_X,&Xwin);
728: PetscLogObjectMemory(draw,sizeof(PetscDraw_X));
729: MPI_Comm_rank(((PetscObject)draw)->comm,&rank);
731: if (!rank) {
732: if (x < 0 || y < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative corner of window");
733: if (w <= 0 || h <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative window width or height");
734: XiQuickWindow(Xwin,draw->display,draw->title,x,y,w,h);
735: MPI_Bcast(&Xwin->win,1,MPI_UNSIGNED_LONG,0,((PetscObject)draw)->comm);
736: } else {
737: unsigned long win = 0;
738: MPI_Bcast(&win,1,MPI_UNSIGNED_LONG,0,((PetscObject)draw)->comm);
739: XiQuickWindowFromWindow(Xwin,draw->display,win);
740: }
742: Xwin->x = x;
743: Xwin->y = y;
744: Xwin->w = w;
745: Xwin->h = h;
746: draw->data = (void*)Xwin;
748: /*
749: Need barrier here so processor 0 does not destroy the window before other
750: processors have completed XiQuickWindow()
751: */
752: PetscDrawClear(draw);
753: PetscDrawSynchronizedFlush(draw);
755: PetscOptionsGetTruth(PETSC_NULL,"-draw_double_buffer",&flg,PETSC_NULL);
756: if (flg) {
757: PetscDrawSetDoubleBuffer(draw);
758: }
760: return(0);
761: }