Actual source code: xops.c

petsc-3.7.7 2017-09-25
Report Typos and Errors
  1: /*
  2:     Defines the operations for the X PetscDraw implementation.
  3: */

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

  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))


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

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

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

 44:   *i = XTRANS(draw,XiWin,x);
 45:   *j = YTRANS(draw,XiWin,y);
 46:   return(0);
 47: }

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

 56:   *x = ITRANS(draw,XiWin,i);
 57:   *y = JTRANS(draw,XiWin,j);
 58:   return(0);
 59: }

 63: static PetscErrorCode PetscDrawPoint_X(PetscDraw draw,PetscReal x,PetscReal y,int c)
 64: {
 65:   PetscDraw_X *XiWin = (PetscDraw_X*)draw->data;
 66:   int         xx,yy,i,j;

 69:   xx = XTRANS(draw,XiWin,x);
 70:   yy = YTRANS(draw,XiWin,y);
 71:   PetscDrawXiSetColor(XiWin,c);
 72:   for (i=-1; i<2; i++) {
 73:     for (j=-1; j<2; j++) {
 74:       XDrawPoint(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,xx+i,yy+j);
 75:     }
 76:   }
 77:   return(0);
 78: }

 82: static PetscErrorCode PetscDrawPointPixel_X(PetscDraw draw,int x,int y,int c)
 83: {
 84:   PetscDraw_X *XiWin = (PetscDraw_X*)draw->data;

 87:   PetscDrawXiSetColor(XiWin,c);
 88:   XDrawPoint(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x,y);
 89:   return(0);
 90: }

 94: static PetscErrorCode PetscDrawLine_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:   return(0);
105: }

109: static PetscErrorCode PetscDrawArrow_X(PetscDraw draw,PetscReal xl,PetscReal yl,PetscReal xr,PetscReal yr,int cl)
110: {
111:   PetscDraw_X *XiWin = (PetscDraw_X*)draw->data;
112:   int         x_1,y_1,x_2,y_2;

115:   PetscDrawXiSetColor(XiWin,cl);
116:   x_1 = XTRANS(draw,XiWin,xl); x_2 = XTRANS(draw,XiWin,xr);
117:   y_1 = YTRANS(draw,XiWin,yl); y_2 = YTRANS(draw,XiWin,yr);
118:   XDrawLine(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x_1,y_1,x_2,y_2);
119:   if (x_1 == x_2 && y_1 == y_2) return(0);
120:   if (x_1 == x_2 && PetscAbs(y_1 - y_2) > 7) {
121:     if (y_2 > y_1) {
122:       XDrawLine(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x_2,y_2,x_2-3,y_2-3);
123:       XDrawLine(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x_2,y_2,x_2+3,y_2-3);
124:     } else {
125:       XDrawLine(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x_2,y_2,x_2-3,y_2+3);
126:       XDrawLine(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x_2,y_2,x_2+3,y_2+3);
127:     }
128:   }
129:   if (y_1 == y_2 && PetscAbs(x_1 - x_2) > 7) {
130:     if (x_2 > x_1) {
131:       XDrawLine(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x_2-3,y_2-3,x_2,y_2);
132:       XDrawLine(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x_2-3,y_2+3,x_2,y_2);
133:     } else {
134:       XDrawLine(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x_2,y_2,x_2+3,y_2-3);
135:       XDrawLine(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x_2,y_2,x_2+3,y_2+3);
136:     }
137:   }
138:   return(0);
139: }

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

149:   PetscDrawXiSetColor(XiWin,c);
150:   x = XTRANS(draw,XiWin,xl); w = XTRANS(draw,XiWin,xr) + 1 - x; if (w <= 0) w = 1;
151:   y = YTRANS(draw,XiWin,yr); h = YTRANS(draw,XiWin,yl) + 1 - y; if (h <= 0) h = 1;
152:   XFillRectangle(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x,y,w,h);
153:   return(0);
154: }

158: static PetscErrorCode PetscDrawEllipse_X(PetscDraw draw,PetscReal x,PetscReal y,PetscReal a,PetscReal b,int c)
159: {
160:   PetscDraw_X *XiWin = (PetscDraw_X*)draw->data;
161:   int         xA,yA,w,h;

164:   PetscDrawXiSetColor(XiWin, c);
165:   xA = XTRANS(draw,XiWin, x - a/2); w = XTRANS(draw,XiWin, x + a/2) + 1 - xA; w = PetscAbs(w);
166:   yA = YTRANS(draw,XiWin, y + b/2); h = YTRANS(draw,XiWin, y - b/2) + 1 - yA; h = PetscAbs(h);
167:   XFillArc(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,xA,yA,w,h,0,360*64);
168:   return(0);
169: }

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

175: 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)
176: {
177:   PetscDraw_X    *XiWin = (PetscDraw_X*)draw->data;

181:   if (c1 == c2 && c2 == c3) {
182:     XPoint pt[3];
183:     PetscDrawXiSetColor(XiWin,c1);
184:     pt[0].x = XTRANS(draw,XiWin,X1);
185:     pt[0].y = YTRANS(draw,XiWin,Y_1);
186:     pt[1].x = XTRANS(draw,XiWin,X2);
187:     pt[1].y = YTRANS(draw,XiWin,Y2);
188:     pt[2].x = XTRANS(draw,XiWin,X3);
189:     pt[2].y = YTRANS(draw,XiWin,Y3);
190:     XFillPolygon(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,pt,3,Convex,CoordModeOrigin);
191:   } else {
192:     int x1,y_1,x2,y2,x3,y3;
193:     x1   = XTRANS(draw,XiWin,X1);
194:     y_1  = YTRANS(draw,XiWin,Y_1);
195:     x2   = XTRANS(draw,XiWin,X2);
196:     y2   = YTRANS(draw,XiWin,Y2);
197:     x3   = XTRANS(draw,XiWin,X3);
198:     y3   = YTRANS(draw,XiWin,Y3);
199:     PetscDrawInterpolatedTriangle_X(XiWin,x1,y_1,c1,x2,y2,c2,x3,y3,c3);
200:   }
201:   return(0);
202: }

206: static PetscErrorCode PetscDrawStringSetSize_X(PetscDraw draw,PetscReal x,PetscReal y)
207: {
208:   PetscDraw_X    *XiWin = (PetscDraw_X*)draw->data;
209:   int            w,h;

213:   w    = (int)((XiWin->w)*x*(draw->port_xr - draw->port_xl)/(draw->coor_xr - draw->coor_xl));
214:   h    = (int)((XiWin->h)*y*(draw->port_yr - draw->port_yl)/(draw->coor_yr - draw->coor_yl));
215:   PetscFree(XiWin->font);
216:   PetscDrawXiFontFixed(XiWin,w,h,&XiWin->font);
217:   return(0);
218: }

222: static PetscErrorCode PetscDrawStringGetSize_X(PetscDraw draw,PetscReal *x,PetscReal  *y)
223: {
224:   PetscDraw_X *XiWin = (PetscDraw_X*)draw->data;
225:   PetscReal   w,h;

228:   w = XiWin->font->font_w; h = XiWin->font->font_h;
229:   if (x) *x = w*(draw->coor_xr - draw->coor_xl)/((XiWin->w)*(draw->port_xr - draw->port_xl));
230:   if (y) *y = h*(draw->coor_yr - draw->coor_yl)/((XiWin->h)*(draw->port_yr - draw->port_yl));
231:   return(0);
232: }

236: static PetscErrorCode PetscDrawString_X(PetscDraw draw,PetscReal x,PetscReal y,int c,const char chrs[])
237: {
238:   PetscDraw_X    *XiWin = (PetscDraw_X*)draw->data;
239:   int            xx,yy,descent = XiWin->font->font_descent;
240:   size_t         len;
241:   char           *substr;
242:   PetscToken     token;

246:   xx = XTRANS(draw,XiWin,x);
247:   yy = YTRANS(draw,XiWin,y);
248:   PetscDrawXiSetColor(XiWin,c);

250:   PetscTokenCreate(chrs,'\n',&token);
251:   PetscTokenFind(token,&substr);
252:   while (substr) {
253:     PetscStrlen(substr,&len);
254:     XDrawString(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,xx,yy-descent,substr,len);
255:     yy  += XiWin->font->font_h;
256:     PetscTokenFind(token,&substr);
257:   }
258:   PetscTokenDestroy(&token);
259:   return(0);
260: }

264: static PetscErrorCode PetscDrawStringVertical_X(PetscDraw draw,PetscReal x,PetscReal y,int c,const char text[])
265: {
266:   PetscDraw_X    *XiWin = (PetscDraw_X*)draw->data;
267:   int            xx,yy,offset = XiWin->font->font_h - XiWin->font->font_descent;
268:   char           chr[2] = {0, 0};

271:   xx = XTRANS(draw,XiWin,x);
272:   yy = YTRANS(draw,XiWin,y);
273:   PetscDrawXiSetColor(XiWin,c);
274:   while ((chr[0] = *text++)) {
275:     XDrawString(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,xx,yy+offset,chr,1);
276:     yy += XiWin->font->font_h;
277:   }
278:   return(0);
279: }

283: static PetscErrorCode PetscDrawFlush_X(PetscDraw draw)
284: {
285:   PetscDraw_X    *XiWin = (PetscDraw_X*)draw->data;
286:   PetscMPIInt    rank;

290:   /* make sure the X server processed requests from all processes */
291:   PetscDrawCollectiveBegin(draw);
292:   XSync(XiWin->disp,False);
293:   PetscDrawCollectiveEnd(draw);
294:   MPI_Barrier(PetscObjectComm((PetscObject)draw));

296:   /* transfer pixmap contents to window (only the first process does this) */
297:   if (XiWin->drw && XiWin->win) {
298:     MPI_Comm_rank(PetscObjectComm((PetscObject)draw),&rank);
299:     PetscDrawCollectiveBegin(draw);
300:     if (!rank) XCopyArea(XiWin->disp,XiWin->drw,XiWin->win,XiWin->gc.set,0,0,XiWin->w,XiWin->h,0,0);
301:     if (!rank) XSync(XiWin->disp,False);
302:     PetscDrawCollectiveEnd(draw);
303:     MPI_Barrier(PetscObjectComm((PetscObject)draw));
304:   }
305:   return(0);
306: }

310: static PetscErrorCode PetscDrawClear_X(PetscDraw draw)
311: {
312:   PetscDraw_X    *XiWin = (PetscDraw_X*)draw->data;
313:   int            xmax = XiWin->w-1,  ymax = XiWin->h-1;
314:   PetscReal      xl = draw->port_xl, yl = draw->port_yl;
315:   PetscReal      xr = draw->port_xr, yr = draw->port_yr;
316:   PetscMPIInt    rank;

320:   /* make sure the X server processed requests from all processes */
321:   PetscDrawCollectiveBegin(draw);
322:   XSync(XiWin->disp,False);
323:   PetscDrawCollectiveEnd(draw);
324:   MPI_Barrier(PetscObjectComm((PetscObject)draw));

326:   /* only the first process handles the clearing business */
327:   PetscDrawCollectiveBegin(draw);
328:   MPI_Comm_rank(PetscObjectComm((PetscObject)draw),&rank);
329:   if (!rank) {
330:     int xa = (int)(xl*xmax), ya = ymax - (int)(yr*ymax);
331:     int xb = (int)(xr*xmax), yb = ymax - (int)(yl*ymax);
332:     unsigned int w = (unsigned int)(xb + 1 - xa);
333:     unsigned int h = (unsigned int)(yb + 1 - ya);
334:     PetscDrawXiSetPixVal(XiWin,XiWin->background);
335:     XFillRectangle(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,xa,ya,w,h);
336:     XSync(XiWin->disp,False);
337:   }
338:   PetscDrawCollectiveEnd(draw);
339:   MPI_Barrier(PetscObjectComm((PetscObject)draw));
340:   return(0);
341: }

345: static PetscErrorCode PetscDrawSetDoubleBuffer_X(PetscDraw draw)
346: {
347:   PetscDraw_X    *win = (PetscDraw_X*)draw->data;
348:   PetscMPIInt    rank;

352:   if (win->drw) return(0);
353:   MPI_Comm_rank(PetscObjectComm((PetscObject)draw),&rank);

355:   PetscDrawCollectiveBegin(draw);
356:   if (!rank) {PetscDrawXiQuickPixmap(win);}
357:   PetscDrawCollectiveEnd(draw);
358:   MPI_Bcast(&win->drw,1,MPI_UNSIGNED_LONG,0,PetscObjectComm((PetscObject)draw));
359:   return(0);
360: }

364: static PetscErrorCode PetscDrawGetPopup_X(PetscDraw draw,PetscDraw *popup)
365: {
366:   PetscDraw_X    *win = (PetscDraw_X*)draw->data;
367:   PetscBool      flg  = PETSC_TRUE;

371:   PetscOptionsGetBool(((PetscObject)draw)->options,((PetscObject)draw)->prefix,"-draw_popup",&flg,NULL);
372:   if (!flg || !win->win) {*popup = NULL; return(0);}

374:   PetscDrawCreate(PetscObjectComm((PetscObject)draw),draw->display,NULL,win->x,win->y+win->h+10,220,220,popup);
375:   PetscObjectSetOptionsPrefix((PetscObject)*popup,"popup_");
376:   PetscObjectAppendOptionsPrefix((PetscObject)*popup,((PetscObject)draw)->prefix);
377:   PetscDrawSetType(*popup,PETSC_DRAW_X);
378:   draw->popup = *popup;
379:   return(0);
380: }

384: static PetscErrorCode PetscDrawSetTitle_X(PetscDraw draw,const char title[])
385: {
386:   PetscDraw_X    *win = (PetscDraw_X*)draw->data;
387:   PetscMPIInt    rank;

391:   if (!win->win) return(0);
392:   MPI_Comm_rank(PetscObjectComm((PetscObject)draw),&rank);
393:   PetscDrawCollectiveBegin(draw);
394:   if (!rank) {
395:     size_t        len;
396:     XTextProperty prop;
397:     PetscStrlen(title,&len);
398:     XGetWMName(win->disp,win->win,&prop);
399:     XFree((void*)prop.value);
400:     prop.value  = (unsigned char*)title;
401:     prop.nitems = (long)len;
402:     XSetWMName(win->disp,win->win,&prop);
403:   }
404:   PetscDrawCollectiveEnd(draw);
405:   return(0);
406: }

410: static PetscErrorCode PetscDrawCheckResizedWindow_X(PetscDraw draw)
411: {
412:   PetscDraw_X    *win = (PetscDraw_X*)draw->data;
413:   int            xywh[4];
414:   PetscMPIInt    rank;

418:   if (!win->win) return(0);
419:   MPI_Comm_rank(PetscObjectComm((PetscObject)draw),&rank);

421:   PetscDrawCollectiveBegin(draw);
422:   if (!rank) {PetscDrawXiGetGeometry(win,xywh,xywh+1,xywh+2,xywh+3);}
423:   PetscDrawCollectiveEnd(draw);
424:   MPI_Bcast(xywh,4,MPI_INT,0,PetscObjectComm((PetscObject)draw));

426:   /* record new window position */
427:   draw->x = win->x = xywh[0];
428:   draw->y = win->y = xywh[1];
429:   if (xywh[2] == win->w && xywh[3] == win->h) return(0);
430:   /* record new window sizes */
431:   draw->w = win->w = xywh[2];
432:   draw->h = win->h = xywh[3];

434:   /* recreate pixmap (only first processor does this) */
435:   PetscDrawCollectiveBegin(draw);
436:   if (!rank && win->drw) {PetscDrawXiQuickPixmap(win);}
437:   PetscDrawCollectiveEnd(draw);
438:   MPI_Bcast(&win->drw,1,MPI_UNSIGNED_LONG,0,PetscObjectComm((PetscObject)draw));
439:   /* reset the clipping */
440:   PetscDrawSetViewport_X(draw,draw->port_xl,draw->port_yl,draw->port_xr,draw->port_yr);
441:   return(0);
442: }

446: static PetscErrorCode PetscDrawResizeWindow_X(PetscDraw draw,int w,int h)
447: {
448:   PetscDraw_X    *win = (PetscDraw_X*)draw->data;
449:   PetscMPIInt    rank;

453:   if (w == win->w && h == win->h) return(0);
454:   MPI_Comm_rank(PetscObjectComm((PetscObject)draw),&rank);

456:   if (win->win) {
457:     PetscDrawCollectiveBegin(draw);
458:     if (!rank) {PetscDrawXiResizeWindow(win,w,h);}
459:     PetscDrawCollectiveEnd(draw);
460:     PetscDrawCheckResizedWindow_X(draw);
461:   } else if (win->drw) {
462:     draw->w = win->w = w; draw->h = win->h = h;
463:     /* recreate pixmap (only first processor does this) */
464:     PetscDrawCollectiveBegin(draw);
465:     if (!rank) {PetscDrawXiQuickPixmap(win);}
466:     MPI_Bcast(&win->drw,1,MPI_UNSIGNED_LONG,0,PetscObjectComm((PetscObject)draw));
467:     /* reset the clipping */
468:     PetscDrawCollectiveEnd(draw);
469:     PetscDrawSetViewport_X(draw,draw->port_xl,draw->port_yl,draw->port_xr,draw->port_yr);
470:   }
471:   return(0);
472: }

474: #include <X11/cursorfont.h>

478: static PetscErrorCode PetscDrawGetMouseButton_X(PetscDraw draw,PetscDrawButton *button,PetscReal *x_user,PetscReal *y_user,PetscReal *x_phys,PetscReal *y_phys)
479: {
480:   PetscDraw_X    *win = (PetscDraw_X*)draw->data;
481:   Cursor         cursor;
482:   XEvent         report;
483:   Window         root,child;
484:   int            root_x,root_y,px=0,py=0;
485:   unsigned int   w,h,border,depth;
486:   unsigned int   keys_button;
487:   PetscMPIInt    rank;
488:   PetscReal      xx,yy;

492:   *button = PETSC_BUTTON_NONE;
493:   if (!win->win) return(0);
494:   MPI_Comm_rank(PetscObjectComm((PetscObject)draw),&rank);

496:   PetscDrawCollectiveBegin(draw);
497:   if (rank) goto finally;

499:   /* change cursor to indicate input */
500:   cursor = XCreateFontCursor(win->disp,XC_hand2); if (!cursor) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Unable to create X cursor");
501:   XDefineCursor(win->disp,win->win,cursor);
502:   /* wait for mouse button events */
503:   XSelectInput(win->disp,win->win,ButtonPressMask|ButtonReleaseMask);
504:   while (XCheckTypedEvent(win->disp,ButtonPress,&report));
505:   XMaskEvent(win->disp,ButtonReleaseMask,&report);
506:   /* get mouse pointer coordinates */
507:   XQueryPointer(win->disp,report.xmotion.window,&root,&child,&root_x,&root_y,&px,&py,&keys_button);
508:   /* the user may resize the window before pressing the mouse button */
509:   XGetGeometry(win->disp,win->win,&root,&root_x,&root_y,&w,&h,&border,&depth);
510:   /* cleanup input event handler and cursor  */
511:   XSelectInput(win->disp,win->win,NoEventMask);
512:   XUndefineCursor(win->disp,win->win);
513:   XFreeCursor(win->disp, cursor);
514:   XSync(win->disp,False);

516:   switch (report.xbutton.button) {
517:   case Button1: *button = PETSC_BUTTON_LEFT; break;
518:   case Button2: *button = PETSC_BUTTON_CENTER; break;
519:   case Button3: *button = PETSC_BUTTON_RIGHT; break;
520:   case Button4: *button = PETSC_BUTTON_WHEEL_UP; break;
521:   case Button5: *button = PETSC_BUTTON_WHEEL_DOWN; break;
522:   }
523:   if (report.xbutton.state & ShiftMask) {
524:     switch (report.xbutton.button) {
525:     case Button1: *button = PETSC_BUTTON_LEFT_SHIFT; break;
526:     case Button2: *button = PETSC_BUTTON_CENTER_SHIFT; break;
527:     case Button3: *button = PETSC_BUTTON_RIGHT_SHIFT; break;
528:     }
529:   }
530:   xx = ((PetscReal)px)/w;
531:   yy = 1 - ((PetscReal)py)/h;
532:   if (x_user) *x_user = draw->coor_xl + (xx - draw->port_xl)*(draw->coor_xr - draw->coor_xl)/(draw->port_xr - draw->port_xl);
533:   if (y_user) *y_user = draw->coor_yl + (yy - draw->port_yl)*(draw->coor_yr - draw->coor_yl)/(draw->port_yr - draw->port_yl);
534:   if (x_phys) *x_phys = xx;
535:   if (y_phys) *y_phys = yy;

537: finally:
538:   PetscDrawCollectiveEnd(draw);
539:   PetscDrawCheckResizedWindow_X(draw);
540:   return(0);
541: }

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

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

563: static PetscErrorCode PetscDrawDestroy_X(PetscDraw draw)
564: {
565:   PetscDraw_X    *win = (PetscDraw_X*)draw->data;

569:   PetscDrawDestroy(&draw->popup);
570:   PetscDrawXiClose(win);
571:   PetscFree(draw->data);
572:   return(0);
573: }

575: static       PetscErrorCode PetscDrawGetSingleton_X(PetscDraw,PetscDraw*);
576: static       PetscErrorCode PetscDrawRestoreSingleton_X(PetscDraw,PetscDraw*);
577: PETSC_INTERN PetscErrorCode PetscDrawGetImage_X(PetscDraw,unsigned char[][3],unsigned int*,unsigned int*,unsigned char*[]);

579: static struct _PetscDrawOps DvOps = { PetscDrawSetDoubleBuffer_X,
580:                                       PetscDrawFlush_X,
581:                                       PetscDrawLine_X,
582:                                       0,
583:                                       0,
584:                                       PetscDrawPoint_X,
585:                                       0,
586:                                       PetscDrawString_X,
587:                                       PetscDrawStringVertical_X,
588:                                       PetscDrawStringSetSize_X,
589:                                       PetscDrawStringGetSize_X,
590:                                       PetscDrawSetViewport_X,
591:                                       PetscDrawClear_X,
592:                                       PetscDrawRectangle_X,
593:                                       PetscDrawTriangle_X,
594:                                       PetscDrawEllipse_X,
595:                                       PetscDrawGetMouseButton_X,
596:                                       PetscDrawPause_X,
597:                                       0,
598:                                       0,
599:                                       PetscDrawGetPopup_X,
600:                                       PetscDrawSetTitle_X,
601:                                       PetscDrawCheckResizedWindow_X,
602:                                       PetscDrawResizeWindow_X,
603:                                       PetscDrawDestroy_X,
604:                                       0,
605:                                       PetscDrawGetSingleton_X,
606:                                       PetscDrawRestoreSingleton_X,
607:                                       0,
608:                                       PetscDrawGetImage_X,
609:                                       0,
610:                                       PetscDrawArrow_X,
611:                                       PetscDrawCoordinateToPixel_X,
612:                                       PetscDrawPixelToCoordinate_X,
613:                                       PetscDrawPointPixel_X,
614:                                       0};


619: static PetscErrorCode PetscDrawGetSingleton_X(PetscDraw draw,PetscDraw *sdraw)
620: {
621:   PetscDraw_X    *Xwin = (PetscDraw_X*)draw->data,*sXwin;

625:   PetscDrawCreate(PETSC_COMM_SELF,draw->display,draw->title,draw->x,draw->y,draw->w,draw->h,sdraw);
626:   PetscObjectChangeTypeName((PetscObject)*sdraw,PETSC_DRAW_X);
627:   PetscMemcpy((*sdraw)->ops,&DvOps,sizeof(DvOps));

629:   if (draw->popup) {
630:     PetscDrawGetSingleton(draw->popup,&(*sdraw)->popup);
631:   }
632:   (*sdraw)->pause   = draw->pause;
633:   (*sdraw)->coor_xl = draw->coor_xl;
634:   (*sdraw)->coor_xr = draw->coor_xr;
635:   (*sdraw)->coor_yl = draw->coor_yl;
636:   (*sdraw)->coor_yr = draw->coor_yr;
637:   (*sdraw)->port_xl = draw->port_xl;
638:   (*sdraw)->port_xr = draw->port_xr;
639:   (*sdraw)->port_yl = draw->port_yl;
640:   (*sdraw)->port_yr = draw->port_yr;

642:   /* share drawables (windows and/or pixmap) from the parent draw */
643:   PetscNewLog(*sdraw,&sXwin);
644:   (*sdraw)->data = (void*)sXwin;
645:   PetscDrawXiInit(sXwin,draw->display);
646:   if (Xwin->win) {
647:     PetscDrawXiQuickWindowFromWindow(sXwin,Xwin->win);
648:     sXwin->drw = Xwin->drw; /* XXX If the window is ever resized, this is wrong! */
649:   } else if (Xwin->drw) {
650:     PetscDrawXiColormap(sXwin);
651:     sXwin->drw = Xwin->drw;
652:   }
653:   PetscDrawXiGetGeometry(sXwin,&sXwin->x,&sXwin->y,&sXwin->w,&sXwin->h);
654:   (*sdraw)->x = sXwin->x; (*sdraw)->y = sXwin->y;
655:   (*sdraw)->w = sXwin->w; (*sdraw)->h = sXwin->h;
656:   return(0);
657: }

661: static PetscErrorCode PetscDrawRestoreSingleton_X(PetscDraw draw,PetscDraw *sdraw)
662: {

666:   if (draw->popup && (*sdraw)->popup) {
667:     PetscBool isdrawx;
668:     PetscDraw_X *pXwin = (PetscDraw_X*)draw->popup->data;
669:     PetscDraw_X *sXwin = (PetscDraw_X*)(*sdraw)->popup->data;
670:     PetscObjectTypeCompare((PetscObject)draw->popup,PETSC_DRAW_X,&isdrawx);
671:     if (!isdrawx) goto finally;
672:     PetscObjectTypeCompare((PetscObject)(*sdraw)->popup,PETSC_DRAW_X,&isdrawx);
673:     if (!isdrawx) goto finally;
674:     if (sXwin->win == pXwin->win) {
675:       PetscDrawRestoreSingleton(draw->popup,&(*sdraw)->popup);
676:     }
677:   }
678: finally:
679:   PetscDrawDestroy(sdraw);
680:   return(0);
681: }

685: static PetscErrorCode PetscDrawXGetDisplaySize_Private(const char name[],int *width,int *height)
686: {
687:   Display *display;

690:   display = XOpenDisplay(name);
691:   if (!display) {
692:     *width  = *height = 0;
693:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Unable to open display on %s\n\
694:     Make sure your COMPUTE NODES are authorized to connect\n\
695:     to this X server and either your DISPLAY variable\n\
696:     is set or you use the -display name option\n",name);
697:   }
698:   *width  = (int)DisplayWidth(display,DefaultScreen(display));
699:   *height = (int)DisplayHeight(display,DefaultScreen(display));
700:   XCloseDisplay(display);
701:   return(0);
702: }

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

707:    Options Database Keys:
708: +  -display <display> - sets the display to use
709: .  -x_virtual - forces use of a X virtual display Xvfb that will not display anything but -draw_save will still work.
710:                 Xvfb is automatically started up in PetscSetDisplay() with this option
711: .  -draw_size w,h - percentage of screeen (either 1, .5, .3, .25), or size in pixels
712: .  -geometry x,y,w,h - set location and size in pixels
713: .  -draw_virtual - do not open a window (draw on a pixmap), -draw_save will still work
714: -  -draw_double_buffer - avoid window flickering (draw on pixmap and flush to window)

716:    Level: beginner

718: .seealso:  PetscDrawOpenX(), PetscDrawSetDisplay(), PetscDrawSetFromOptions()

720: M*/

724: PETSC_EXTERN PetscErrorCode PetscDrawCreate_X(PetscDraw draw)
725: {
726:   PetscDraw_X    *Xwin;
728:   PetscMPIInt    rank;
729:   int            x = draw->x,y = draw->y,w = draw->w,h = draw->h;
730:   static int     xavailable = 0,yavailable = 0,ybottom = 0,xmax = 0,ymax = 0;
731:   PetscBool      set,dvirtual = PETSC_FALSE,doublebuffer = PETSC_TRUE;
732:   PetscInt       xywh[4],osize = 4,nsizes=2;
733:   PetscReal      sizes[2] = {.3,.3};

736:   /* get the display variable */
737:   if (!draw->display) {
738:     PetscMalloc1(256,&draw->display);
739:     PetscGetDisplay(draw->display,256);
740:   }

742:   /* initialize the display size */
743:   if (!xmax) {
744:     PetscDrawXGetDisplaySize_Private(draw->display,&xmax,&ymax);
745:     /* if some processors fail on this and others succed then this is a problem ! */
746:     if (ierr) {
747:       (*PetscErrorPrintf)("PETSc unable to use X windows\nproceeding without graphics\n");
748:       PetscDrawSetType(draw,PETSC_DRAW_NULL);
749:       return(0);
750:     }
751:   }

753:   /* allow user to set size of drawable */
754:   PetscOptionsGetRealArray(((PetscObject)draw)->options,((PetscObject)draw)->prefix,"-draw_size",sizes,&nsizes,&set);
755:   if (set && nsizes == 1 && sizes[0] > 1.0) sizes[1] = sizes[0];
756:   if (set) {
757:     if (sizes[0] > 1.0)       w = (int)sizes[0];
758:     else if (sizes[0] == 1.0) w = PETSC_DRAW_FULL_SIZE;
759:     else if (sizes[0] == .5)  w = PETSC_DRAW_HALF_SIZE;
760:     else if (sizes[0] == .3)  w = PETSC_DRAW_THIRD_SIZE;
761:     else if (sizes[0] == .25) w = PETSC_DRAW_QUARTER_SIZE;
762:     if (sizes[1] > 1.0)       h = (int)sizes[1];
763:     else if (sizes[1] == 1.0) h = PETSC_DRAW_FULL_SIZE;
764:     else if (sizes[1] == .5)  h = PETSC_DRAW_HALF_SIZE;
765:     else if (sizes[1] == .3)  h = PETSC_DRAW_THIRD_SIZE;
766:     else if (sizes[1] == .25) h = PETSC_DRAW_QUARTER_SIZE;
767:   }
768:   if (w == PETSC_DECIDE || w == PETSC_DEFAULT) w = draw->w = 300;
769:   if (h == PETSC_DECIDE || h == PETSC_DEFAULT) h = draw->h = 300;
770:   switch (w) {
771:   case PETSC_DRAW_FULL_SIZE:    w = draw->w = (xmax - 10);   break;
772:   case PETSC_DRAW_HALF_SIZE:    w = draw->w = (xmax - 20)/2; break;
773:   case PETSC_DRAW_THIRD_SIZE:   w = draw->w = (xmax - 30)/3; break;
774:   case PETSC_DRAW_QUARTER_SIZE: w = draw->w = (xmax - 40)/4; break;
775:   }
776:   switch (h) {
777:   case PETSC_DRAW_FULL_SIZE:    h = draw->h = (ymax - 10);   break;
778:   case PETSC_DRAW_HALF_SIZE:    h = draw->h = (ymax - 20)/2; break;
779:   case PETSC_DRAW_THIRD_SIZE:   h = draw->h = (ymax - 30)/3; break;
780:   case PETSC_DRAW_QUARTER_SIZE: h = draw->h = (ymax - 40)/4; break;
781:   }

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

785:   if (!dvirtual) {

787:     /* allow user to set location and size of window */
788:     xywh[0] = x; xywh[1] = y; xywh[2] = w; xywh[3] = h;
789:     PetscOptionsGetIntArray(((PetscObject)draw)->options,((PetscObject)draw)->prefix,"-geometry",xywh,&osize,NULL);
790:     x = (int)xywh[0]; y = (int)xywh[1]; w = (int)xywh[2]; h = (int)xywh[3];
791:     if (w == PETSC_DECIDE || w == PETSC_DEFAULT) w = 300;
792:     if (h == PETSC_DECIDE || h == PETSC_DEFAULT) h = 300;
793:     draw->x = x; draw->y = y; draw->w = w; draw->h = h;

795:     if (draw->x == PETSC_DECIDE || draw->y == PETSC_DECIDE) {
796:       /*
797:        PETSc tries to place windows starting in the upper left corner
798:         and moving across to the right.

800:        +0,0-------------------------------------------+
801:        |  Region used so far  +xavailable,yavailable  |
802:        |                      |                       |
803:        |                      |                       |
804:        +--------------------- +ybottom                |
805:        |                                              |
806:        |                                              |
807:        +----------------------------------------------+xmax,ymax

809:       */
810:       /*  First: can we add it to the right? */
811:       if (xavailable + w + 10 <= xmax) {
812:         x       = xavailable;
813:         y       = yavailable;
814:         ybottom = PetscMax(ybottom,y + h + 30);
815:       } else {
816:         /* No, so add it below on the left */
817:         xavailable = x = 0;
818:         yavailable = y = ybottom;
819:         ybottom    = ybottom + h + 30;
820:       }
821:     }
822:     /* update available region */
823:     xavailable = PetscMax(xavailable,x + w + 10);
824:     if (xavailable >= xmax) {
825:       xavailable = 0;
826:       yavailable = yavailable + h + 30;
827:       ybottom    = yavailable;
828:     }
829:     if (yavailable >= ymax) {
830:       y          = 0;
831:       yavailable = 0;
832:       ybottom    = 0;
833:     }

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

837:   MPI_Comm_rank(PetscObjectComm((PetscObject)draw),&rank);
838:   if (!rank && (w <= 0 || h <= 0)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative window width or height");

840:   PetscNewLog(draw,&Xwin);
841:   PetscMemcpy(draw->ops,&DvOps,sizeof(DvOps));
842:   draw->data = (void*)Xwin;

844:   PetscDrawXiInit(Xwin,draw->display);
845:   if (!dvirtual) {
846:     Xwin->x = x; Xwin->y = y;
847:     Xwin->w = w; Xwin->h = h;
848:     if (!rank) {PetscDrawXiQuickWindow(Xwin,draw->title,x,y,w,h);}
849:     MPI_Bcast(&Xwin->win,1,MPI_UNSIGNED_LONG,0,PetscObjectComm((PetscObject)draw));
850:     if (rank) {PetscDrawXiQuickWindowFromWindow(Xwin,Xwin->win);}
851:   } else {
852:     Xwin->x = 0; Xwin->y = 0;
853:     Xwin->w = w; Xwin->h = h;
854:     PetscDrawXiColormap(Xwin);
855:     if (!rank) {PetscDrawXiQuickPixmap(Xwin);}
856:     MPI_Bcast(&Xwin->drw,1,MPI_UNSIGNED_LONG,0,PetscObjectComm((PetscObject)draw));
857:   }
858:   PetscDrawXiGetGeometry(Xwin,&Xwin->x,&Xwin->y,&Xwin->w,&Xwin->h);
859:   draw->x = Xwin->x; draw->y = Xwin->y;
860:   draw->w = Xwin->w; draw->h = Xwin->h;

862:   PetscOptionsGetBool(((PetscObject)draw)->options,((PetscObject)draw)->prefix,"-draw_double_buffer",&doublebuffer,NULL);
863:   if (doublebuffer) {PetscDrawSetDoubleBuffer(draw);}
864:   return(0);
865: }