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;

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

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

 37:   *i = XTRANS(draw,XiWin,x);
 38:   *j = YTRANS(draw,XiWin,y);
 39:   return 0;
 40: }

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

 46:   *x = ITRANS(draw,XiWin,i);
 47:   *y = JTRANS(draw,XiWin,j);
 48:   return 0;
 49: }

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

 56:   xx = XTRANS(draw,XiWin,x);
 57:   yy = YTRANS(draw,XiWin,y);
 58:   PetscDrawXiSetColor(XiWin,c);
 59:   for (i=-1; i<2; i++) {
 60:     for (j=-1; j<2; j++) {
 61:       XDrawPoint(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,xx+i,yy+j);
 62:     }
 63:   }
 64:   return 0;
 65: }

 67: static PetscErrorCode PetscDrawPointPixel_X(PetscDraw draw,int x,int y,int c)
 68: {
 69:   PetscDraw_X *XiWin = (PetscDraw_X*)draw->data;

 71:   PetscDrawXiSetColor(XiWin,c);
 72:   XDrawPoint(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x,y);
 73:   return 0;
 74: }

 76: static PetscErrorCode PetscDrawLine_X(PetscDraw draw,PetscReal xl,PetscReal yl,PetscReal xr,PetscReal yr,int cl)
 77: {
 78:   PetscDraw_X *XiWin = (PetscDraw_X*)draw->data;
 79:   int         x_1,y_1,x_2,y_2;

 81:   PetscDrawXiSetColor(XiWin,cl);
 82:   x_1 = XTRANS(draw,XiWin,xl); x_2  = XTRANS(draw,XiWin,xr);
 83:   y_1 = YTRANS(draw,XiWin,yl); y_2  = YTRANS(draw,XiWin,yr);
 84:   XDrawLine(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x_1,y_1,x_2,y_2);
 85:   return 0;
 86: }

 88: static PetscErrorCode PetscDrawArrow_X(PetscDraw draw,PetscReal xl,PetscReal yl,PetscReal xr,PetscReal yr,int cl)
 89: {
 90:   PetscDraw_X *XiWin = (PetscDraw_X*)draw->data;
 91:   int         x_1,y_1,x_2,y_2;

 93:   PetscDrawXiSetColor(XiWin,cl);
 94:   x_1 = XTRANS(draw,XiWin,xl); x_2 = XTRANS(draw,XiWin,xr);
 95:   y_1 = YTRANS(draw,XiWin,yl); y_2 = YTRANS(draw,XiWin,yr);
 96:   XDrawLine(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x_1,y_1,x_2,y_2);
 97:   if (x_1 == x_2 && y_1 == y_2) return 0;
 98:   if (x_1 == x_2 && PetscAbs(y_1 - y_2) > 7) {
 99:     if (y_2 > y_1) {
100:       XDrawLine(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x_2,y_2,x_2-3,y_2-3);
101:       XDrawLine(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x_2,y_2,x_2+3,y_2-3);
102:     } else {
103:       XDrawLine(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x_2,y_2,x_2-3,y_2+3);
104:       XDrawLine(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x_2,y_2,x_2+3,y_2+3);
105:     }
106:   }
107:   if (y_1 == y_2 && PetscAbs(x_1 - x_2) > 7) {
108:     if (x_2 > x_1) {
109:       XDrawLine(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x_2-3,y_2-3,x_2,y_2);
110:       XDrawLine(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x_2-3,y_2+3,x_2,y_2);
111:     } else {
112:       XDrawLine(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x_2,y_2,x_2+3,y_2-3);
113:       XDrawLine(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x_2,y_2,x_2+3,y_2+3);
114:     }
115:   }
116:   return 0;
117: }

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

124:   PetscDrawXiSetColor(XiWin,c);
125:   x = XTRANS(draw,XiWin,xl); w = XTRANS(draw,XiWin,xr) + 1 - x; if (w <= 0) w = 1;
126:   y = YTRANS(draw,XiWin,yr); h = YTRANS(draw,XiWin,yl) + 1 - y; if (h <= 0) h = 1;
127:   XFillRectangle(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x,y,w,h);
128:   return 0;
129: }

131: static PetscErrorCode PetscDrawEllipse_X(PetscDraw draw,PetscReal x,PetscReal y,PetscReal a,PetscReal b,int c)
132: {
133:   PetscDraw_X *XiWin = (PetscDraw_X*)draw->data;
134:   int         xA,yA,w,h;

136:   PetscDrawXiSetColor(XiWin, c);
137:   xA = XTRANS(draw,XiWin, x - a/2); w = XTRANS(draw,XiWin, x + a/2) + 1 - xA; w = PetscAbs(w);
138:   yA = YTRANS(draw,XiWin, y + b/2); h = YTRANS(draw,XiWin, y - b/2) + 1 - yA; h = PetscAbs(h);
139:   XFillArc(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,xA,yA,w,h,0,360*64);
140:   return 0;
141: }

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

145: 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)
146: {
147:   PetscDraw_X *XiWin = (PetscDraw_X*)draw->data;

149:   if (c1 == c2 && c2 == c3) {
150:     XPoint pt[3];
151:     PetscDrawXiSetColor(XiWin,c1);
152:     pt[0].x = XTRANS(draw,XiWin,X1);
153:     pt[0].y = YTRANS(draw,XiWin,Y_1);
154:     pt[1].x = XTRANS(draw,XiWin,X2);
155:     pt[1].y = YTRANS(draw,XiWin,Y2);
156:     pt[2].x = XTRANS(draw,XiWin,X3);
157:     pt[2].y = YTRANS(draw,XiWin,Y3);
158:     XFillPolygon(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,pt,3,Convex,CoordModeOrigin);
159:   } else {
160:     int x1,y_1,x2,y2,x3,y3;
161:     x1   = XTRANS(draw,XiWin,X1);
162:     y_1  = YTRANS(draw,XiWin,Y_1);
163:     x2   = XTRANS(draw,XiWin,X2);
164:     y2   = YTRANS(draw,XiWin,Y2);
165:     x3   = XTRANS(draw,XiWin,X3);
166:     y3   = YTRANS(draw,XiWin,Y3);
167:     PetscDrawInterpolatedTriangle_X(XiWin,x1,y_1,c1,x2,y2,c2,x3,y3,c3);
168:   }
169:   return 0;
170: }

172: static PetscErrorCode PetscDrawStringSetSize_X(PetscDraw draw,PetscReal x,PetscReal y)
173: {
174:   PetscDraw_X *XiWin = (PetscDraw_X*)draw->data;
175:   int          w,h;

177:   w    = (int)((XiWin->w)*x*(draw->port_xr - draw->port_xl)/(draw->coor_xr - draw->coor_xl));
178:   h    = (int)((XiWin->h)*y*(draw->port_yr - draw->port_yl)/(draw->coor_yr - draw->coor_yl));
179:   PetscFree(XiWin->font);
180:   PetscDrawXiFontFixed(XiWin,w,h,&XiWin->font);
181:   return 0;
182: }

184: static PetscErrorCode PetscDrawStringGetSize_X(PetscDraw draw,PetscReal *x,PetscReal  *y)
185: {
186:   PetscDraw_X *XiWin = (PetscDraw_X*)draw->data;
187:   PetscReal   w,h;

189:   w = XiWin->font->font_w; h = XiWin->font->font_h;
190:   if (x) *x = w*(draw->coor_xr - draw->coor_xl)/((XiWin->w)*(draw->port_xr - draw->port_xl));
191:   if (y) *y = h*(draw->coor_yr - draw->coor_yl)/((XiWin->h)*(draw->port_yr - draw->port_yl));
192:   return 0;
193: }

195: static PetscErrorCode PetscDrawString_X(PetscDraw draw,PetscReal x,PetscReal y,int c,const char chrs[])
196: {
197:   PetscDraw_X *XiWin         = (PetscDraw_X*)draw->data;
198:   int          xx,yy,descent = XiWin->font->font_descent;
199:   size_t       len;
200:   char        *substr;
201:   PetscToken   token;

203:   xx = XTRANS(draw,XiWin,x);
204:   yy = YTRANS(draw,XiWin,y);
205:   PetscDrawXiSetColor(XiWin,c);

207:   PetscTokenCreate(chrs,'\n',&token);
208:   PetscTokenFind(token,&substr);
209:   while (substr) {
210:     PetscStrlen(substr,&len);
211:     XDrawString(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,xx,yy-descent,substr,len);
212:     yy  += XiWin->font->font_h;
213:     PetscTokenFind(token,&substr);
214:   }
215:   PetscTokenDestroy(&token);
216:   return 0;
217: }

219: static PetscErrorCode PetscDrawStringVertical_X(PetscDraw draw,PetscReal x,PetscReal y,int c,const char text[])
220: {
221:   PetscDraw_X    *XiWin = (PetscDraw_X*)draw->data;
222:   int            xx,yy,offset = XiWin->font->font_h - XiWin->font->font_descent;
223:   char           chr[2] = {0, 0};

225:   xx = XTRANS(draw,XiWin,x);
226:   yy = YTRANS(draw,XiWin,y);
227:   PetscDrawXiSetColor(XiWin,c);
228:   while ((chr[0] = *text++)) {
229:     XDrawString(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,xx,yy+offset,chr,1);
230:     yy += XiWin->font->font_h;
231:   }
232:   return 0;
233: }

235: static PetscErrorCode PetscDrawFlush_X(PetscDraw draw)
236: {
237:   PetscDraw_X    *XiWin = (PetscDraw_X*)draw->data;
238:   PetscMPIInt    rank;

241:   /* make sure the X server processed requests from all processes */
242:   PetscDrawCollectiveBegin(draw);
243:   XSync(XiWin->disp,False);
244:   PetscDrawCollectiveEnd(draw);
245:   MPI_Barrier(PetscObjectComm((PetscObject)draw));

247:   /* transfer pixmap contents to window (only the first process does this) */
248:   if (XiWin->drw && XiWin->win) {
249:     MPI_Comm_rank(PetscObjectComm((PetscObject)draw),&rank);
250:     PetscDrawCollectiveBegin(draw);
251:     if (rank == 0) XCopyArea(XiWin->disp,XiWin->drw,XiWin->win,XiWin->gc.set,0,0,XiWin->w,XiWin->h,0,0);
252:     if (rank == 0) XSync(XiWin->disp,False);
253:     PetscDrawCollectiveEnd(draw);
254:     MPI_Barrier(PetscObjectComm((PetscObject)draw));
255:   }
256:   return 0;
257: }

259: static PetscErrorCode PetscDrawClear_X(PetscDraw draw)
260: {
261:   PetscDraw_X    *XiWin = (PetscDraw_X*)draw->data;
262:   int            xmax = XiWin->w-1,  ymax = XiWin->h-1;
263:   PetscReal      xl = draw->port_xl, yl = draw->port_yl;
264:   PetscReal      xr = draw->port_xr, yr = draw->port_yr;
265:   PetscMPIInt    rank;

268:   /* make sure the X server processed requests from all processes */
269:   PetscDrawCollectiveBegin(draw);
270:   XSync(XiWin->disp,False);
271:   PetscDrawCollectiveEnd(draw);
272:   MPI_Barrier(PetscObjectComm((PetscObject)draw));

274:   /* only the first process handles the clearing business */
275:   PetscDrawCollectiveBegin(draw);
276:   MPI_Comm_rank(PetscObjectComm((PetscObject)draw),&rank);
277:   if (rank == 0) {
278:     int xa = (int)(xl*xmax), ya = ymax - (int)(yr*ymax);
279:     int xb = (int)(xr*xmax), yb = ymax - (int)(yl*ymax);
280:     unsigned int w = (unsigned int)(xb + 1 - xa);
281:     unsigned int h = (unsigned int)(yb + 1 - ya);
282:     PetscDrawXiSetPixVal(XiWin,XiWin->background);
283:     XFillRectangle(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,xa,ya,w,h);
284:     XSync(XiWin->disp,False);
285:   }
286:   PetscDrawCollectiveEnd(draw);
287:   MPI_Barrier(PetscObjectComm((PetscObject)draw));
288:   return 0;
289: }

291: static PetscErrorCode PetscDrawSetDoubleBuffer_X(PetscDraw draw)
292: {
293:   PetscDraw_X    *win = (PetscDraw_X*)draw->data;
294:   PetscMPIInt    rank;

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

300:   PetscDrawCollectiveBegin(draw);
301:   if (rank == 0) PetscDrawXiQuickPixmap(win);
302:   PetscDrawCollectiveEnd(draw);
303:   MPI_Bcast(&win->drw,1,MPI_UNSIGNED_LONG,0,PetscObjectComm((PetscObject)draw));
304:   return 0;
305: }

307: static PetscErrorCode PetscDrawGetPopup_X(PetscDraw draw,PetscDraw *popup)
308: {
309:   PetscDraw_X *win = (PetscDraw_X*)draw->data;
310:   PetscBool    flg = PETSC_TRUE;

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

315:   PetscDrawCreate(PetscObjectComm((PetscObject)draw),draw->display,NULL,win->x,win->y+win->h+10,220,220,popup);
316:   PetscObjectSetOptionsPrefix((PetscObject)*popup,"popup_");
317:   PetscObjectAppendOptionsPrefix((PetscObject)*popup,((PetscObject)draw)->prefix);
318:   PetscDrawSetType(*popup,PETSC_DRAW_X);
319:   draw->popup = *popup;
320:   return 0;
321: }

323: static PetscErrorCode PetscDrawSetTitle_X(PetscDraw draw,const char title[])
324: {
325:   PetscDraw_X    *win = (PetscDraw_X*)draw->data;
326:   PetscMPIInt    rank;

329:   if (!win->win) return 0;
330:   MPI_Comm_rank(PetscObjectComm((PetscObject)draw),&rank);
331:   PetscDrawCollectiveBegin(draw);
332:   if (rank == 0) {
333:     size_t        len;
334:     XTextProperty prop;
335:     PetscStrlen(title,&len);
336:     XGetWMName(win->disp,win->win,&prop);
337:     XFree((void*)prop.value);
338:     prop.value  = (unsigned char*)title;
339:     prop.nitems = (long)len;
340:     XSetWMName(win->disp,win->win,&prop);
341:   }
342:   PetscDrawCollectiveEnd(draw);
343:   return 0;
344: }

346: static PetscErrorCode PetscDrawCheckResizedWindow_X(PetscDraw draw)
347: {
348:   PetscDraw_X    *win = (PetscDraw_X*)draw->data;
349:   int            xywh[4];
350:   PetscMPIInt    rank;

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

356:   PetscDrawCollectiveBegin(draw);
357:   if (rank == 0) PetscDrawXiGetGeometry(win,xywh,xywh+1,xywh+2,xywh+3);
358:   PetscDrawCollectiveEnd(draw);
359:   MPI_Bcast(xywh,4,MPI_INT,0,PetscObjectComm((PetscObject)draw));

361:   /* record new window position */
362:   draw->x = win->x = xywh[0];
363:   draw->y = win->y = xywh[1];
364:   if (xywh[2] == win->w && xywh[3] == win->h) return 0;
365:   /* record new window sizes */
366:   draw->w = win->w = xywh[2];
367:   draw->h = win->h = xywh[3];

369:   /* recreate pixmap (only first processor does this) */
370:   PetscDrawCollectiveBegin(draw);
371:   if (rank == 0 && win->drw) PetscDrawXiQuickPixmap(win);
372:   PetscDrawCollectiveEnd(draw);
373:   MPI_Bcast(&win->drw,1,MPI_UNSIGNED_LONG,0,PetscObjectComm((PetscObject)draw));
374:   /* reset the clipping */
375:   PetscDrawSetViewport_X(draw,draw->port_xl,draw->port_yl,draw->port_xr,draw->port_yr);
376:   return 0;
377: }

379: static PetscErrorCode PetscDrawResizeWindow_X(PetscDraw draw,int w,int h)
380: {
381:   PetscDraw_X    *win = (PetscDraw_X*)draw->data;
382:   PetscMPIInt    rank;

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

388:   if (win->win) {
389:     PetscDrawCollectiveBegin(draw);
390:     if (rank == 0) PetscDrawXiResizeWindow(win,w,h);
391:     PetscDrawCollectiveEnd(draw);
392:     PetscDrawCheckResizedWindow_X(draw);
393:   } else if (win->drw) {
394:     draw->w = win->w = w; draw->h = win->h = h;
395:     /* recreate pixmap (only first processor does this) */
396:     PetscDrawCollectiveBegin(draw);
397:     if (rank == 0) PetscDrawXiQuickPixmap(win);
398:     MPI_Bcast(&win->drw,1,MPI_UNSIGNED_LONG,0,PetscObjectComm((PetscObject)draw));
399:     /* reset the clipping */
400:     PetscDrawCollectiveEnd(draw);
401:     PetscDrawSetViewport_X(draw,draw->port_xl,draw->port_yl,draw->port_xr,draw->port_yr);
402:   }
403:   return 0;
404: }

406: #include <X11/cursorfont.h>

408: static PetscErrorCode PetscDrawGetMouseButton_X(PetscDraw draw,PetscDrawButton *button,PetscReal *x_user,PetscReal *y_user,PetscReal *x_phys,PetscReal *y_phys)
409: {
410:   PetscDraw_X    *win = (PetscDraw_X*)draw->data;
411:   Cursor         cursor;
412:   XEvent         report;
413:   Window         root,child;
414:   int            root_x,root_y,px=0,py=0;
415:   unsigned int   w,h,border,depth;
416:   unsigned int   keys_button;
417:   PetscMPIInt    rank;
418:   PetscReal      xx,yy;

421:   *button = PETSC_BUTTON_NONE;
422:   if (!win->win) return 0;
423:   MPI_Comm_rank(PetscObjectComm((PetscObject)draw),&rank);

425:   PetscDrawCollectiveBegin(draw);
426:   if (rank) goto finally;

428:   /* change cursor to indicate input */
430:   XDefineCursor(win->disp,win->win,cursor);
431:   /* wait for mouse button events */
432:   XSelectInput(win->disp,win->win,ButtonPressMask|ButtonReleaseMask);
433:   while (XCheckTypedEvent(win->disp,ButtonPress,&report));
434:   XMaskEvent(win->disp,ButtonReleaseMask,&report);
435:   /* get mouse pointer coordinates */
436:   XQueryPointer(win->disp,report.xmotion.window,&root,&child,&root_x,&root_y,&px,&py,&keys_button);
437:   /* the user may resize the window before pressing the mouse button */
438:   XGetGeometry(win->disp,win->win,&root,&root_x,&root_y,&w,&h,&border,&depth);
439:   /* cleanup input event handler and cursor  */
440:   XSelectInput(win->disp,win->win,NoEventMask);
441:   XUndefineCursor(win->disp,win->win);
442:   XFreeCursor(win->disp, cursor);
443:   XSync(win->disp,False);

445:   switch (report.xbutton.button) {
446:   case Button1: *button = PETSC_BUTTON_LEFT; break;
447:   case Button2: *button = PETSC_BUTTON_CENTER; break;
448:   case Button3: *button = PETSC_BUTTON_RIGHT; break;
449:   case Button4: *button = PETSC_BUTTON_WHEEL_UP; break;
450:   case Button5: *button = PETSC_BUTTON_WHEEL_DOWN; break;
451:   }
452:   if (report.xbutton.state & ShiftMask) {
453:     switch (report.xbutton.button) {
454:     case Button1: *button = PETSC_BUTTON_LEFT_SHIFT; break;
455:     case Button2: *button = PETSC_BUTTON_CENTER_SHIFT; break;
456:     case Button3: *button = PETSC_BUTTON_RIGHT_SHIFT; break;
457:     }
458:   }
459:   xx = ((PetscReal)px)/w;
460:   yy = 1 - ((PetscReal)py)/h;
461:   if (x_user) *x_user = draw->coor_xl + (xx - draw->port_xl)*(draw->coor_xr - draw->coor_xl)/(draw->port_xr - draw->port_xl);
462:   if (y_user) *y_user = draw->coor_yl + (yy - draw->port_yl)*(draw->coor_yr - draw->coor_yl)/(draw->port_yr - draw->port_yl);
463:   if (x_phys) *x_phys = xx;
464:   if (y_phys) *y_phys = yy;

466: finally:
467:   PetscDrawCollectiveEnd(draw);
468:   PetscDrawCheckResizedWindow_X(draw);
469:   return 0;
470: }

472: static PetscErrorCode PetscDrawPause_X(PetscDraw draw)
473: {
474:   PetscDraw_X *win = (PetscDraw_X*)draw->data;

476:   if (!win->win) return 0;
477:   if (draw->pause > 0) PetscSleep(draw->pause);
478:   else if (draw->pause == -1) {
479:     PetscDrawButton button = PETSC_BUTTON_NONE;
480:     PetscDrawGetMouseButton(draw,&button,NULL,NULL,NULL,NULL);
481:     if (button == PETSC_BUTTON_CENTER) draw->pause = 0;
482:   }
483:   return 0;
484: }

486: static PetscErrorCode PetscDrawDestroy_X(PetscDraw draw)
487: {
488:   PetscDraw_X *win = (PetscDraw_X*)draw->data;

490:   PetscDrawDestroy(&draw->popup);
491:   PetscDrawXiClose(win);
492:   PetscFree(draw->data);
493:   return 0;
494: }

496: static       PetscErrorCode PetscDrawGetSingleton_X(PetscDraw,PetscDraw*);
497: static       PetscErrorCode PetscDrawRestoreSingleton_X(PetscDraw,PetscDraw*);
498: PETSC_INTERN PetscErrorCode PetscDrawGetImage_X(PetscDraw,unsigned char[][3],unsigned int*,unsigned int*,unsigned char*[]);

500: static struct _PetscDrawOps DvOps = { PetscDrawSetDoubleBuffer_X,
501:                                       PetscDrawFlush_X,
502:                                       PetscDrawLine_X,
503:                                       NULL,
504:                                       NULL,
505:                                       PetscDrawPoint_X,
506:                                       NULL,
507:                                       PetscDrawString_X,
508:                                       PetscDrawStringVertical_X,
509:                                       PetscDrawStringSetSize_X,
510:                                       PetscDrawStringGetSize_X,
511:                                       PetscDrawSetViewport_X,
512:                                       PetscDrawClear_X,
513:                                       PetscDrawRectangle_X,
514:                                       PetscDrawTriangle_X,
515:                                       PetscDrawEllipse_X,
516:                                       PetscDrawGetMouseButton_X,
517:                                       PetscDrawPause_X,
518:                                       NULL,
519:                                       NULL,
520:                                       PetscDrawGetPopup_X,
521:                                       PetscDrawSetTitle_X,
522:                                       PetscDrawCheckResizedWindow_X,
523:                                       PetscDrawResizeWindow_X,
524:                                       PetscDrawDestroy_X,
525:                                       NULL,
526:                                       PetscDrawGetSingleton_X,
527:                                       PetscDrawRestoreSingleton_X,
528:                                       NULL,
529:                                       PetscDrawGetImage_X,
530:                                       NULL,
531:                                       PetscDrawArrow_X,
532:                                       PetscDrawCoordinateToPixel_X,
533:                                       PetscDrawPixelToCoordinate_X,
534:                                       PetscDrawPointPixel_X,
535:                                       NULL};

537: static PetscErrorCode PetscDrawGetSingleton_X(PetscDraw draw,PetscDraw *sdraw)
538: {
539:   PetscDraw_X *Xwin = (PetscDraw_X*)draw->data,*sXwin;

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

545:   if (draw->popup) PetscDrawGetSingleton(draw->popup,&(*sdraw)->popup);
546:   (*sdraw)->pause   = draw->pause;
547:   (*sdraw)->coor_xl = draw->coor_xl;
548:   (*sdraw)->coor_xr = draw->coor_xr;
549:   (*sdraw)->coor_yl = draw->coor_yl;
550:   (*sdraw)->coor_yr = draw->coor_yr;
551:   (*sdraw)->port_xl = draw->port_xl;
552:   (*sdraw)->port_xr = draw->port_xr;
553:   (*sdraw)->port_yl = draw->port_yl;
554:   (*sdraw)->port_yr = draw->port_yr;

556:   /* share drawables (windows and/or pixmap) from the parent draw */
557:   PetscNewLog(*sdraw,&sXwin);
558:   (*sdraw)->data = (void*)sXwin;
559:   PetscDrawXiInit(sXwin,draw->display);
560:   if (Xwin->win) {
561:     PetscDrawXiQuickWindowFromWindow(sXwin,Xwin->win);
562:     sXwin->drw = Xwin->drw; /* XXX If the window is ever resized, this is wrong! */
563:   } else if (Xwin->drw) {
564:     PetscDrawXiColormap(sXwin);
565:     sXwin->drw = Xwin->drw;
566:   }
567:   PetscDrawXiGetGeometry(sXwin,&sXwin->x,&sXwin->y,&sXwin->w,&sXwin->h);
568:   (*sdraw)->x = sXwin->x; (*sdraw)->y = sXwin->y;
569:   (*sdraw)->w = sXwin->w; (*sdraw)->h = sXwin->h;
570:   return 0;
571: }

573: static PetscErrorCode PetscDrawRestoreSingleton_X(PetscDraw draw,PetscDraw *sdraw)
574: {
575:   if (draw->popup && (*sdraw)->popup) {
576:     PetscBool isdrawx;
577:     PetscDraw_X *pXwin = (PetscDraw_X*)draw->popup->data;
578:     PetscDraw_X *sXwin = (PetscDraw_X*)(*sdraw)->popup->data;
579:     PetscObjectTypeCompare((PetscObject)draw->popup,PETSC_DRAW_X,&isdrawx);
580:     if (!isdrawx) goto finally;
581:     PetscObjectTypeCompare((PetscObject)(*sdraw)->popup,PETSC_DRAW_X,&isdrawx);
582:     if (!isdrawx) goto finally;
583:     if (sXwin->win == pXwin->win) {
584:       PetscDrawRestoreSingleton(draw->popup,&(*sdraw)->popup);
585:     }
586:   }
587: finally:
588:   PetscDrawDestroy(sdraw);
589:   return 0;
590: }

592: static PetscErrorCode PetscDrawXGetDisplaySize_Private(const char name[],int *width,int *height,PetscBool *has_display)
593: {
594:   Display *display;

596:   display = XOpenDisplay(name);
597:   if (!display) {
598:     *width  = *height = 0;
599:     (*PetscErrorPrintf)("Unable to open display on %s\n\
600:     Make sure your COMPUTE NODES are authorized to connect\n\
601:     to this X server and either your DISPLAY variable\n\
602:     is set or you use the -display name option\n",name);
603:     *has_display = PETSC_FALSE;
604:     return 0;
605:   }
606:   *has_display = PETSC_TRUE;
607:   *width  = (int)DisplayWidth(display,DefaultScreen(display));
608:   *height = (int)DisplayHeight(display,DefaultScreen(display));
609:   XCloseDisplay(display);
610:   return 0;
611: }

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

616:    Options Database Keys:
617: +  -display <display> - sets the display to use
618: .  -x_virtual - forces use of a X virtual display Xvfb that will not display anything but -draw_save will still work.
619:                 Xvfb is automatically started up in PetscSetDisplay() with this option
620: .  -draw_size w,h - percentage of screen (either 1, .5, .3, .25), or size in pixels
621: .  -geometry x,y,w,h - set location and size in pixels
622: .  -draw_virtual - do not open a window (draw on a pixmap), -draw_save will still work
623: -  -draw_double_buffer - avoid window flickering (draw on pixmap and flush to window)

625:    Level: beginner

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

629: M*/

631: PETSC_EXTERN PetscErrorCode PetscDrawCreate_X(PetscDraw draw)
632: {
633:   PetscDraw_X    *Xwin;
634:   PetscMPIInt    rank;
635:   int            x = draw->x,y = draw->y,w = draw->w,h = draw->h;
636:   static int     xavailable = 0,yavailable = 0,ybottom = 0,xmax = 0,ymax = 0;
637:   PetscBool      set,dvirtual = PETSC_FALSE,doublebuffer = PETSC_TRUE,has_display;
638:   PetscInt       xywh[4],osize = 4,nsizes=2;
639:   PetscReal      sizes[2] = {.3,.3};
640:   static size_t  DISPLAY_LENGTH = 265;

642:   /* get the display variable */
643:   if (!draw->display) {
644:     PetscMalloc1(DISPLAY_LENGTH,&draw->display);
645:     PetscGetDisplay(draw->display,DISPLAY_LENGTH);
646:   }

648:   /* initialize the display size */
649:   if (!xmax) {
650:     PetscDrawXGetDisplaySize_Private(draw->display,&xmax,&ymax,&has_display);
651:     /* if some processors fail on this and others succed then this is a problem ! */
652:     if (!has_display) {
653:       (*PetscErrorPrintf)("PETSc unable to use X windows\nproceeding without graphics\n");
654:       PetscDrawSetType(draw,PETSC_DRAW_NULL);
655:       return 0;
656:     }
657:   }

659:   /* allow user to set size of drawable */
660:   PetscOptionsGetRealArray(((PetscObject)draw)->options,((PetscObject)draw)->prefix,"-draw_size",sizes,&nsizes,&set);
661:   if (set && nsizes == 1 && sizes[0] > 1.0) sizes[1] = sizes[0];
662:   if (set) {
663:     if (sizes[0] > 1.0)       w = (int)sizes[0];
664:     else if (sizes[0] == 1.0) w = PETSC_DRAW_FULL_SIZE;
665:     else if (sizes[0] == .5)  w = PETSC_DRAW_HALF_SIZE;
666:     else if (sizes[0] == .3)  w = PETSC_DRAW_THIRD_SIZE;
667:     else if (sizes[0] == .25) w = PETSC_DRAW_QUARTER_SIZE;
668:     if (sizes[1] > 1.0)       h = (int)sizes[1];
669:     else if (sizes[1] == 1.0) h = PETSC_DRAW_FULL_SIZE;
670:     else if (sizes[1] == .5)  h = PETSC_DRAW_HALF_SIZE;
671:     else if (sizes[1] == .3)  h = PETSC_DRAW_THIRD_SIZE;
672:     else if (sizes[1] == .25) h = PETSC_DRAW_QUARTER_SIZE;
673:   }
674:   if (w == PETSC_DECIDE || w == PETSC_DEFAULT) w = draw->w = 300;
675:   if (h == PETSC_DECIDE || h == PETSC_DEFAULT) h = draw->h = 300;
676:   switch (w) {
677:   case PETSC_DRAW_FULL_SIZE:    w = draw->w = (xmax - 10);   break;
678:   case PETSC_DRAW_HALF_SIZE:    w = draw->w = (xmax - 20)/2; break;
679:   case PETSC_DRAW_THIRD_SIZE:   w = draw->w = (xmax - 30)/3; break;
680:   case PETSC_DRAW_QUARTER_SIZE: w = draw->w = (xmax - 40)/4; break;
681:   }
682:   switch (h) {
683:   case PETSC_DRAW_FULL_SIZE:    h = draw->h = (ymax - 10);   break;
684:   case PETSC_DRAW_HALF_SIZE:    h = draw->h = (ymax - 20)/2; break;
685:   case PETSC_DRAW_THIRD_SIZE:   h = draw->h = (ymax - 30)/3; break;
686:   case PETSC_DRAW_QUARTER_SIZE: h = draw->h = (ymax - 40)/4; break;
687:   }

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

691:   if (!dvirtual) {

693:     /* allow user to set location and size of window */
694:     xywh[0] = x; xywh[1] = y; xywh[2] = w; xywh[3] = h;
695:     PetscOptionsGetIntArray(((PetscObject)draw)->options,((PetscObject)draw)->prefix,"-geometry",xywh,&osize,NULL);
696:     x = (int)xywh[0]; y = (int)xywh[1]; w = (int)xywh[2]; h = (int)xywh[3];
697:     if (w == PETSC_DECIDE || w == PETSC_DEFAULT) w = 300;
698:     if (h == PETSC_DECIDE || h == PETSC_DEFAULT) h = 300;
699:     draw->x = x; draw->y = y; draw->w = w; draw->h = h;

701:     if (draw->x == PETSC_DECIDE || draw->y == PETSC_DECIDE) {
702:       /*
703:        PETSc tries to place windows starting in the upper left corner
704:         and moving across to the right.

706:        +0,0-------------------------------------------+
707:        |  Region used so far  +xavailable,yavailable  |
708:        |                      |                       |
709:        |                      |                       |
710:        +--------------------- +ybottom                |
711:        |                                              |
712:        |                                              |
713:        +----------------------------------------------+xmax,ymax

715:       */
716:       /*  First: can we add it to the right? */
717:       if (xavailable + w + 10 <= xmax) {
718:         x       = xavailable;
719:         y       = yavailable;
720:         ybottom = PetscMax(ybottom,y + h + 30);
721:       } else {
722:         /* No, so add it below on the left */
723:         xavailable = x = 0;
724:         yavailable = y = ybottom;
725:         ybottom    = ybottom + h + 30;
726:       }
727:     }
728:     /* update available region */
729:     xavailable = PetscMax(xavailable,x + w + 10);
730:     if (xavailable >= xmax) {
731:       xavailable = 0;
732:       yavailable = yavailable + h + 30;
733:       ybottom    = yavailable;
734:     }
735:     if (yavailable >= ymax) {
736:       y          = 0;
737:       yavailable = 0;
738:       ybottom    = 0;
739:     }

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

743:   MPI_Comm_rank(PetscObjectComm((PetscObject)draw),&rank);

746:   PetscNewLog(draw,&Xwin);
747:   PetscMemcpy(draw->ops,&DvOps,sizeof(DvOps));
748:   draw->data = (void*)Xwin;

750:   PetscDrawXiInit(Xwin,draw->display);
751:   if (!dvirtual) {
752:     Xwin->x = x; Xwin->y = y;
753:     Xwin->w = w; Xwin->h = h;
754:     if (rank == 0) PetscDrawXiQuickWindow(Xwin,draw->title,x,y,w,h);
755:     MPI_Bcast(&Xwin->win,1,MPI_UNSIGNED_LONG,0,PetscObjectComm((PetscObject)draw));
756:     if (rank) PetscDrawXiQuickWindowFromWindow(Xwin,Xwin->win);
757:   } else {
758:     Xwin->x = 0; Xwin->y = 0;
759:     Xwin->w = w; Xwin->h = h;
760:     PetscDrawXiColormap(Xwin);
761:     if (rank == 0) PetscDrawXiQuickPixmap(Xwin);
762:     MPI_Bcast(&Xwin->drw,1,MPI_UNSIGNED_LONG,0,PetscObjectComm((PetscObject)draw));
763:   }
764:   PetscDrawXiGetGeometry(Xwin,&Xwin->x,&Xwin->y,&Xwin->w,&Xwin->h);
765:   draw->x = Xwin->x; draw->y = Xwin->y;
766:   draw->w = Xwin->w; draw->h = Xwin->h;

768:   PetscOptionsGetBool(((PetscObject)draw)->options,((PetscObject)draw)->prefix,"-draw_double_buffer",&doublebuffer,NULL);
769:   if (doublebuffer) PetscDrawSetDoubleBuffer(draw);
770:   return 0;
771: }