Actual source code: xops.c

petsc-3.6.4 2016-04-12
Report Typos and Errors
  2: /*
  3:     Defines the operations for the X PetscDraw implementation.
  4: */

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

  8: /*
  9:      These macros transform from the users coordinates to the  X-window pixel coordinates.
 10: */
 11: #define XTRANS(draw,xwin,x)  (int)(((xwin)->w)*((draw)->port_xl + (((x - (draw)->coor_xl)*((draw)->port_xr - (draw)->port_xl))/((draw)->coor_xr - (draw)->coor_xl))))
 12: #define YTRANS(draw,xwin,y)  (int)(((xwin)->h)*(1.0-(draw)->port_yl - (((y - (draw)->coor_yl)*((draw)->port_yr - (draw)->port_yl))/((draw)->coor_yr - (draw)->coor_yl))))

 14: #define ITRANS(draw,xwin,i)  (draw)->coor_xl + (i*((draw)->coor_xr - (draw)->coor_xl)/((xwin)->w) - (draw)->port_xl)/((draw)->port_xr - (draw)->port_xl)
 15: #define JTRANS(draw,xwin,j)  draw->coor_yl + (((double)j)/xwin->h + draw->port_yl - 1.0)*(draw->coor_yr - draw->coor_yl)/(draw->port_yl - draw->port_yr)

 19: PetscErrorCode PetscDrawCoordinateToPixel_X(PetscDraw draw,PetscReal x,PetscReal y,PetscInt *i,PetscInt *j)
 20: {
 21:   PetscDraw_X *XiWin = (PetscDraw_X*)draw->data;

 24:   *i = XTRANS(draw,XiWin,x);
 25:   *j = YTRANS(draw,XiWin,y);
 26:   return(0);
 27: }

 31: PetscErrorCode PetscDrawPixelToCoordinate_X(PetscDraw draw,PetscInt i,PetscInt j,PetscReal *x,PetscReal *y)
 32: {
 33:   PetscDraw_X *XiWin = (PetscDraw_X*)draw->data;

 36:   *x = ITRANS(draw,XiWin,i);
 37:   *y = JTRANS(draw,XiWin,j);
 38:   return(0);
 39: }

 43: PetscErrorCode PetscDrawLine_X(PetscDraw draw,PetscReal xl,PetscReal yl,PetscReal xr,PetscReal yr,int cl)
 44: {
 45:   PetscDraw_X *XiWin = (PetscDraw_X*)draw->data;
 46:   int         x1,y_1,x2,y2;

 49:   PetscDrawXiSetColor(XiWin,cl);
 50:   x1  = XTRANS(draw,XiWin,xl);   x2  = XTRANS(draw,XiWin,xr);
 51:   y_1 = YTRANS(draw,XiWin,yl);   y2  = YTRANS(draw,XiWin,yr);
 52:   if (x1 == x2 && y_1 == y2) return(0);
 53:   XDrawLine(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x1,y_1,x2,y2);
 54:   return(0);
 55: }

 59: PetscErrorCode PetscDrawArrow_X(PetscDraw draw,PetscReal xl,PetscReal yl,PetscReal xr,PetscReal yr,int cl)
 60: {
 61:   PetscDraw_X *XiWin = (PetscDraw_X*)draw->data;
 62:   int         x1,y_1,x2,y2;

 65:   PetscDrawXiSetColor(XiWin,cl);
 66:   x1  = XTRANS(draw,XiWin,xl);   x2  = XTRANS(draw,XiWin,xr);
 67:   y_1 = YTRANS(draw,XiWin,yl);   y2  = YTRANS(draw,XiWin,yr);
 68:   if (x1 == x2 && y_1 == y2) return(0);
 69:   XDrawLine(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x1,y_1,x2,y2);
 70:   if (x1 == x2 && PetscAbs(y_1 - y2) > 7) {
 71:     if (y2 > y_1) {
 72:       XDrawLine(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x2,y2,x2-3,y2-3);
 73:       XDrawLine(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x2,y2,x2+3,y2-3);
 74:     } else {
 75:       XDrawLine(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x2,y2,x2-3,y2+3);
 76:       XDrawLine(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x2,y2,x2+3,y2+3);
 77:     }
 78:   }
 79:   return(0);
 80: }

 84: static PetscErrorCode PetscDrawPoint_X(PetscDraw draw,PetscReal x,PetscReal y,int c)
 85: {
 86:   int         xx,yy,i,j;
 87:   PetscDraw_X *XiWin = (PetscDraw_X*)draw->data;

 90:   xx = XTRANS(draw,XiWin,x);  yy = YTRANS(draw,XiWin,y);
 91:   PetscDrawXiSetColor(XiWin,c);
 92:   for (i=-1; i<2; i++) {
 93:     for (j=-1; j<2; j++) {
 94:       XDrawPoint(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,xx+i,yy+j);
 95:     }
 96:   }
 97:   return(0);
 98: }

102: static PetscErrorCode PetscDrawPointPixel_X(PetscDraw draw,PetscInt x,PetscInt y,int c)
103: {
104:   PetscDraw_X *XiWin = (PetscDraw_X*)draw->data;

107:   PetscDrawXiSetColor(XiWin,c);
108:   XDrawPoint(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x,y);
109:   return(0);
110: }

114: static PetscErrorCode PetscDrawRectangle_X(PetscDraw draw,PetscReal xl,PetscReal yl,PetscReal xr,PetscReal yr,int c1,int c2,int c3,int c4)
115: {
116:   PetscDraw_X *XiWin = (PetscDraw_X*)draw->data;
117:   int         x1,y_1,w,h,c = (c1 + c2 + c3 + c4)/4;

120:   PetscDrawXiSetColor(XiWin,c);
121:   x1  = XTRANS(draw,XiWin,xl);   w  = XTRANS(draw,XiWin,xr) - x1;
122:   y_1 = YTRANS(draw,XiWin,yr);   h  = YTRANS(draw,XiWin,yl) - y_1;
123:   if (w <= 0) w = 1;
124:   if (h <= 0) h = 1;
125:   XFillRectangle(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x1,y_1,w,h);
126:   return(0);
127: }

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

137:   PetscDrawXiSetColor(XiWin, c);
138:   xA = XTRANS(Win, XiWin, x - a/2.0); w = XTRANS(Win, XiWin, x + a/2.0) - xA;
139:   yA = YTRANS(Win, XiWin, y + b/2.0); h = PetscAbs(YTRANS(Win, XiWin, y - b/2.0) - yA);
140:   XFillArc(XiWin->disp, PetscDrawXiDrawable(XiWin), XiWin->gc.set, xA, yA, w, h, 0, 23040);
141:   return(0);
142: }

144: extern PetscErrorCode PetscDrawInterpolatedTriangle_X(PetscDraw_X*,int,int,int,int,int,int,int,int,int);

148: 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)
149: {
150:   PetscDraw_X    *XiWin = (PetscDraw_X*)draw->data;

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

179: static PetscErrorCode PetscDrawString_X(PetscDraw draw,PetscReal x,PetscReal y,int c,const char chrs[])
180: {
182:   int            xx,yy;
183:   size_t         len;
184:   PetscDraw_X    *XiWin = (PetscDraw_X*)draw->data;
185:   char           *substr;
186:   PetscToken     token;

189:   xx = XTRANS(draw,XiWin,x);  yy = YTRANS(draw,XiWin,y);
190:   PetscDrawXiSetColor(XiWin,c);

192:   PetscTokenCreate(chrs,'\n',&token);
193:   PetscTokenFind(token,&substr);
194:   PetscStrlen(substr,&len);
195:   XDrawString(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,xx,yy - XiWin->font->font_descent,substr,len);
196:   PetscTokenFind(token,&substr);
197:   while (substr) {
198:     yy  += 4*XiWin->font->font_descent;
199:     PetscStrlen(substr,&len);
200:     XDrawString(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,xx,yy - XiWin->font->font_descent,substr,len);
201:     PetscTokenFind(token,&substr);
202:   }
203:   PetscTokenDestroy(&token);
204:   return(0);
205: }

207: extern PetscErrorCode PetscDrawXiFontFixed(PetscDraw_X*,int,int,PetscDrawXiFont**);

211: static PetscErrorCode PetscDrawStringSetSize_X(PetscDraw draw,PetscReal x,PetscReal y)
212: {
213:   PetscDraw_X    *XiWin = (PetscDraw_X*)draw->data;
215:   int            w,h;

218:   w    = (int)((XiWin->w)*x*(draw->port_xr - draw->port_xl)/(draw->coor_xr - draw->coor_xl));
219:   h    = (int)((XiWin->h)*y*(draw->port_yr - draw->port_yl)/(draw->coor_yr - draw->coor_yl));
220:   PetscFree(XiWin->font);
221:   PetscDrawXiFontFixed(XiWin,w,h,&XiWin->font);
222:   return(0);
223: }

227: PetscErrorCode PetscDrawStringGetSize_X(PetscDraw draw,PetscReal *x,PetscReal  *y)
228: {
229:   PetscDraw_X *XiWin = (PetscDraw_X*)draw->data;
230:   PetscReal   w,h;

233:   w = XiWin->font->font_w; h = XiWin->font->font_h;
234:   if (x) *x = w*(draw->coor_xr - draw->coor_xl)/((XiWin->w)*(draw->port_xr - draw->port_xl));
235:   if (y) *y = h*(draw->coor_yr - draw->coor_yl)/((XiWin->h)*(draw->port_yr - draw->port_yl));
236:   return(0);
237: }

241: PetscErrorCode PetscDrawStringVertical_X(PetscDraw draw,PetscReal x,PetscReal y,int c,const char chrs[])
242: {
244:   int            xx,yy;
245:   PetscDraw_X    *XiWin = (PetscDraw_X*)draw->data;
246:   char           tmp[2];
247:   PetscReal      tw,th;
248:   size_t         i,n;

251:   PetscStrlen(chrs,&n);
252:   tmp[1] = 0;
253:   PetscDrawXiSetColor(XiWin,c);
254:   PetscDrawStringGetSize_X(draw,&tw,&th);
255:   xx   = XTRANS(draw,XiWin,x);
256:   for (i=0; i<n; i++) {
257:     tmp[0] = chrs[i];
258:     yy     = YTRANS(draw,XiWin,y-th*i);
259:     XDrawString(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set, xx,yy - XiWin->font->font_descent,tmp,1);
260:   }
261:   return(0);
262: }

266: static PetscErrorCode PetscDrawFlush_X(PetscDraw draw)
267: {
268:   PetscDraw_X*   XiWin = (PetscDraw_X*)draw->data;

272:   if (XiWin->drw && XiWin->win) XCopyArea(XiWin->disp,XiWin->drw,XiWin->win,XiWin->gc.set,0,0,XiWin->w,XiWin->h,0,0);
273:   XFlush(XiWin->disp);
274:   XSync(XiWin->disp,False);
275:   if (draw->saveonflush) {PetscDrawSave(draw);}
276:   return(0);
277: }

281: static PetscErrorCode PetscDrawSynchronizedFlush_X(PetscDraw draw)
282: {
284:   PetscMPIInt    rank;
285:   PetscDraw_X    *XiWin = (PetscDraw_X*)draw->data;

288:   if (XiWin->drw && XiWin->win) {
289:     MPI_Comm_rank(PetscObjectComm((PetscObject)draw),&rank);
290:     /* make sure data has actually arrived at server */
291:     XSync(XiWin->disp,False);
292:     MPI_Barrier(PetscObjectComm((PetscObject)draw));
293:     if (!rank) {
294:       XCopyArea(XiWin->disp,XiWin->drw,XiWin->win,XiWin->gc.set,0,0,XiWin->w,XiWin->h,0,0);
295:       XFlush(XiWin->disp);
296:       XSync(XiWin->disp,False);
297:     }
298:     MPI_Barrier(PetscObjectComm((PetscObject)draw));
299:   } else {
300:     MPI_Barrier(PetscObjectComm((PetscObject)draw));
301:     XSync(XiWin->disp,False);
302:     MPI_Barrier(PetscObjectComm((PetscObject)draw));
303:   }
304:   return(0);
305: }

309: static PetscErrorCode PetscDrawSetViewport_X(PetscDraw draw,PetscReal xl,PetscReal yl,PetscReal xr,PetscReal yr)
310: {
311:   PetscDraw_X *XiWin = (PetscDraw_X*)draw->data;
312:   XRectangle  box;

315:   box.x     = (int)(xl*XiWin->w);     box.y      = (int)((1.0-yr)*XiWin->h);
316:   box.width = (int)((xr-xl)*XiWin->w);box.height = (int)((yr-yl)*XiWin->h);
317:   XSetClipRectangles(XiWin->disp,XiWin->gc.set,0,0,&box,1,Unsorted);
318:   return(0);
319: }

323: static PetscErrorCode PetscDrawClear_X(PetscDraw draw)
324: {
325:   PetscDraw_X    *XiWin = (PetscDraw_X*)draw->data;
326:   int            x, y, w, h;

330:   PetscDrawSave(draw);
331:   x    = (int)(draw->port_xl*XiWin->w);
332:   w    = (int)((draw->port_xr - draw->port_xl)*XiWin->w);
333:   y    = (int)((1.0-draw->port_yr)*XiWin->h);
334:   h    = (int)((draw->port_yr - draw->port_yl)*XiWin->h);
335:   PetscDrawXiSetPixVal(XiWin,XiWin->background);
336:   XFillRectangle(XiWin->disp,PetscDrawXiDrawable(XiWin),XiWin->gc.set,x,y,w,h);
337:   return(0);
338: }

342: static PetscErrorCode PetscDrawSynchronizedClear_X(PetscDraw draw)
343: {
345:   PetscMPIInt    rank;
346:   PetscDraw_X    *XiWin = (PetscDraw_X*)draw->data;

349:   MPI_Barrier(PetscObjectComm((PetscObject)draw));
350:   MPI_Comm_rank(PetscObjectComm((PetscObject)draw),&rank);
351:   if (!rank) {
352:     PetscDrawClear_X(draw);
353:   }
354:   XFlush(XiWin->disp);
355:   XFlush(XiWin->disp);
356:   MPI_Barrier(PetscObjectComm((PetscObject)draw));
357:   /*  XSync(XiWin->disp,False); */
358:   MPI_Barrier(PetscObjectComm((PetscObject)draw));
359:   return(0);
360: }

364: static PetscErrorCode PetscDrawSetDoubleBuffer_X(PetscDraw draw)
365: {
366:   PetscDraw_X    *win = (PetscDraw_X*)draw->data;
368:   PetscMPIInt    rank;

371:   if (win->drw) return(0);

373:   MPI_Comm_rank(PetscObjectComm((PetscObject)draw),&rank);
374:   if (!rank) win->drw = XCreatePixmap(win->disp,win->win,win->w,win->h,win->depth);

376:   /* try to make sure it is actually done before passing info to all */
377:   XSync(win->disp,False);
378:   MPI_Bcast(&win->drw,1,MPI_UNSIGNED_LONG,0,PetscObjectComm((PetscObject)draw));
379:   return(0);
380: }

382: #include <X11/cursorfont.h>

386: static PetscErrorCode PetscDrawGetMouseButton_X(PetscDraw draw,PetscDrawButton *button,PetscReal *x_user,PetscReal *y_user,PetscReal *x_phys,PetscReal *y_phys)
387: {
388:   XEvent       report;
389:   PetscDraw_X  *win = (PetscDraw_X*)draw->data;
390:   Window       root,child;
391:   int          root_x,root_y,px,py;
392:   unsigned int keys_button;
393:   Cursor       cursor = 0;

396:   /* change cursor to indicate input */
397:   if (!cursor) {
398:     cursor = XCreateFontCursor(win->disp,XC_hand2);
399:     if (!cursor) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Unable to create X cursor");
400:   }
401:   XDefineCursor(win->disp,win->win,cursor);
402:   XSelectInput(win->disp,win->win,ButtonPressMask | ButtonReleaseMask);

404:   while (XCheckTypedEvent(win->disp,ButtonPress,&report));
405:   XMaskEvent(win->disp,ButtonReleaseMask,&report);
406:   switch (report.xbutton.button) {
407:   case Button1:
408:     if (report.xbutton.state & ShiftMask) *button = PETSC_BUTTON_LEFT_SHIFT;
409:     else                                  *button = PETSC_BUTTON_LEFT;
410:     break;
411:   case Button2:
412:     if (report.xbutton.state & ShiftMask) *button = PETSC_BUTTON_CENTER_SHIFT;
413:     else                                  *button = PETSC_BUTTON_CENTER;
414:     break;
415:   case Button3:
416:     if (report.xbutton.state & ShiftMask) *button = PETSC_BUTTON_RIGHT_SHIFT;
417:     else                                  *button = PETSC_BUTTON_RIGHT;
418:     break;
419:   }
420:   XQueryPointer(win->disp,report.xmotion.window,&root,&child,&root_x,&root_y,&px,&py,&keys_button);

422:   if (x_phys) *x_phys = ((double)px)/((double)win->w);
423:   if (y_phys) *y_phys = 1.0 - ((double)py)/((double)win->h);

425:   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);
426:   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);

428:   XUndefineCursor(win->disp,win->win);
429:   XFlush(win->disp);
430:   XSync(win->disp,False);
431:   return(0);
432: }

436: static PetscErrorCode PetscDrawPause_X(PetscDraw draw)
437: {
439:   PetscDraw_X    *win = (PetscDraw_X*)draw->data;

442:   if (!win->win) return(0);
443:   if (draw->pause > 0) PetscSleep(draw->pause);
444:   else if (draw->pause == -1) {
445:     PetscDrawButton button;
446:     PetscMPIInt     rank;
447:     MPI_Comm_rank(PetscObjectComm((PetscObject)draw),&rank);
448:     if (!rank) {
449:       PetscDrawGetMouseButton(draw,&button,0,0,0,0);
450:       if (button == PETSC_BUTTON_CENTER) draw->pause = 0;
451:     }
452:     MPI_Bcast(&draw->pause,1,MPIU_REAL,0,PetscObjectComm((PetscObject)draw));
453:   }
454:   return(0);
455: }

459: static PetscErrorCode PetscDrawGetPopup_X(PetscDraw draw,PetscDraw *popup)
460: {
462:   PetscDraw_X    *win = (PetscDraw_X*)draw->data;
463:   PetscBool      flg   = PETSC_TRUE;

466:   PetscOptionsGetBool(((PetscObject)draw)->prefix,"-draw_popup",&flg,NULL);
467:   if (flg) {
468:     PetscDrawOpenX(PetscObjectComm((PetscObject)draw),NULL,NULL,win->x,win->y+win->h+36,220,220,popup);
469:     draw->popup = *popup;
470:   } else {
471:     *popup = NULL;
472:   }
473:   return(0);
474: }

478: static PetscErrorCode PetscDrawSetTitle_X(PetscDraw draw,const char title[])
479: {
480:   PetscDraw_X    *win = (PetscDraw_X*)draw->data;
481:   XTextProperty  prop;
483:   size_t         len;

486:   if (win->win) {
487:     XGetWMName(win->disp,win->win,&prop);
488:     XFree((void*)prop.value);
489:     prop.value  = (unsigned char*)title;
490:     PetscStrlen(title,&len);
491:     prop.nitems = (long) len;
492:     XSetWMName(win->disp,win->win,&prop);
493:   }
494:   return(0);
495: }

499: static PetscErrorCode PetscDrawResizeWindow_X(PetscDraw draw,int w,int h)
500: {
501:   PetscDraw_X    *win = (PetscDraw_X*)draw->data;
502:   unsigned int   ww,hh,border,depth;
503:   int            x,y;
505:   Window         root;

508:   if (win->win) {
509:     XResizeWindow(win->disp,win->win,w,h);
510:     XGetGeometry(win->disp,win->win,&root,&x,&y,&ww,&hh,&border,&depth);
511:     PetscDrawCheckResizedWindow(draw);
512:   }
513:   return(0);
514: }

518: static PetscErrorCode PetscDrawCheckResizedWindow_X(PetscDraw draw)
519: {
520:   PetscDraw_X    *win = (PetscDraw_X*)draw->data;
522:   int            x,y;
523:   PetscMPIInt    rank;
524:   Window         root;
525:   unsigned int   w,h,border,depth,geo[2];
526:   PetscReal      xl,xr,yl,yr;
527:   XRectangle     box;

530:   if (!win->win) return(0);
531:   MPI_Comm_rank(PetscObjectComm((PetscObject)draw),&rank);
532:   if (!rank) {
533:     XFlush(win->disp);
534:     XSync(win->disp,False);
535:     XSync(win->disp,False);
536:     XGetGeometry(win->disp,win->win,&root,&x,&y,geo,geo+1,&border,&depth);
537:   }
538:   MPI_Bcast(geo,2,MPI_UNSIGNED,0,PetscObjectComm((PetscObject)draw));
539:   w    = geo[0];
540:   h    = geo[1];
541:   if (w == (unsigned int) win->w && h == (unsigned int) win->h) return(0);

543:   /* record new window sizes */

545:   win->h = h; win->w = w;

547:   /* Free buffer space and create new version (only first processor does this) */
548:   if (win->drw) win->drw = XCreatePixmap(win->disp,win->win,win->w,win->h,win->depth);

550:   /* reset the clipping */
551:   xl = draw->port_xl; yl = draw->port_yl;
552:   xr = draw->port_xr; yr = draw->port_yr;

554:   box.x     = (int)(xl*win->w);     box.y      = (int)((1.0-yr)*win->h);
555:   box.width = (int)((xr-xl)*win->w);box.height = (int)((yr-yl)*win->h);
556:   XSetClipRectangles(win->disp,win->gc.set,0,0,&box,1,Unsorted);

558:   /* try to make sure it is actually done before passing info to all */
559:   XSync(win->disp,False);
560:   MPI_Bcast(&win->drw,1,MPI_UNSIGNED_LONG,0,PetscObjectComm((PetscObject)draw));
561:   return(0);
562: }

564: static PetscErrorCode PetscDrawGetSingleton_X(PetscDraw,PetscDraw*);
565: static PetscErrorCode PetscDrawRestoreSingleton_X(PetscDraw,PetscDraw*);

569: PetscErrorCode PetscDrawDestroy_X(PetscDraw draw)
570: {
571:   PetscDraw_X    *win = (PetscDraw_X*)draw->data;
573: #if defined(PETSC_HAVE_POPEN)
574:   char           command[PETSC_MAX_PATH_LEN];
575:   PetscMPIInt    rank;
576:   FILE           *fd;
577: #endif

580:   if (draw->savefinalfilename) {
581:     PetscDrawSetSave(draw,draw->savefinalfilename,PETSC_FALSE);
582:     draw->savefilecount = 0;
583:   }
584:   PetscDrawSynchronizedClear(draw);

586: #if defined(PETSC_HAVE_POPEN)
587:   MPI_Comm_rank(PetscObjectComm((PetscObject)draw),&rank);
588:   if (draw->savefilename && !rank && draw->savefilemovie) {
589:     PetscSNPrintf(command,PETSC_MAX_PATH_LEN,"ffmpeg  -i %s/%s_%%d.Gif %s.m4v",draw->savefilename,draw->savefilename,draw->savefilename);
590:     PetscPOpen(PETSC_COMM_SELF,NULL,command,"r",&fd);
591:     PetscPClose(PETSC_COMM_SELF,fd,NULL);
592:   }
593: #endif
594:   PetscBarrier((PetscObject)draw);

596:   XFreeGC(win->disp,win->gc.set);
597:   XCloseDisplay(win->disp);
598:   PetscDrawDestroy(&draw->popup);
599:   PetscFree(win->font);
600:   PetscFree(draw->data);
601:   return(0);
602: }

604: PetscErrorCode PetscDrawSave_X(PetscDraw);
605: PetscErrorCode PetscDrawSetSave_X(PetscDraw,const char*);

607: static struct _PetscDrawOps DvOps = { PetscDrawSetDoubleBuffer_X,
608:                                       PetscDrawFlush_X,
609:                                       PetscDrawLine_X,
610:                                       0,
611:                                       0,
612:                                       PetscDrawPoint_X,
613:                                       0,
614:                                       PetscDrawString_X,
615:                                       PetscDrawStringVertical_X,
616:                                       PetscDrawStringSetSize_X,
617:                                       PetscDrawStringGetSize_X,
618:                                       PetscDrawSetViewport_X,
619:                                       PetscDrawClear_X,
620:                                       PetscDrawSynchronizedFlush_X,
621:                                       PetscDrawRectangle_X,
622:                                       PetscDrawTriangle_X,
623:                                       PetscDrawEllipse_X,
624:                                       PetscDrawGetMouseButton_X,
625:                                       PetscDrawPause_X,
626:                                       PetscDrawSynchronizedClear_X,
627:                                       0,
628:                                       0,
629:                                       PetscDrawGetPopup_X,
630:                                       PetscDrawSetTitle_X,
631:                                       PetscDrawCheckResizedWindow_X,
632:                                       PetscDrawResizeWindow_X,
633:                                       PetscDrawDestroy_X,
634:                                       0,
635:                                       PetscDrawGetSingleton_X,
636:                                       PetscDrawRestoreSingleton_X,
637: #if defined(PETSC_HAVE_AFTERIMAGE)
638:                                       PetscDrawSave_X,
639: #else
640:                                       0,
641: #endif
642:                                       PetscDrawSetSave_X,
643:                                       0,
644:                                       PetscDrawArrow_X,
645:                                       PetscDrawCoordinateToPixel_X,
646:                                       PetscDrawPixelToCoordinate_X,
647:                                       PetscDrawPointPixel_X,
648:                                       0};


651: extern PetscErrorCode PetscDrawXiQuickWindow(PetscDraw_X*,char*,char*,int,int,int,int);
652: extern PetscErrorCode PetscDrawXiQuickWindowFromWindow(PetscDraw_X*,char*,Window);

656: static PetscErrorCode PetscDrawGetSingleton_X(PetscDraw draw,PetscDraw *sdraw)
657: {
659:   PetscDraw_X    *Xwin = (PetscDraw_X*)draw->data,*sXwin;

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

666:   (*sdraw)->ops->destroy = NULL;

668:   (*sdraw)->pause   = draw->pause;
669:   (*sdraw)->coor_xl = draw->coor_xl;
670:   (*sdraw)->coor_xr = draw->coor_xr;
671:   (*sdraw)->coor_yl = draw->coor_yl;
672:   (*sdraw)->coor_yr = draw->coor_yr;
673:   (*sdraw)->port_xl = draw->port_xl;
674:   (*sdraw)->port_xr = draw->port_xr;
675:   (*sdraw)->port_yl = draw->port_yl;
676:   (*sdraw)->port_yr = draw->port_yr;
677:   (*sdraw)->popup   = draw->popup;

679:   /* actually create and open the window */
680:   PetscNew(&sXwin);
681:   PetscDrawXiQuickWindowFromWindow(sXwin,draw->display,Xwin->win);

683:   sXwin->x       = Xwin->x;
684:   sXwin->y       = Xwin->y;
685:   sXwin->w       = Xwin->w;
686:   sXwin->h       = Xwin->h;
687:   (*sdraw)->data = (void*)sXwin;
688:   return(0);
689: }

693: static PetscErrorCode PetscDrawRestoreSingleton_X(PetscDraw draw,PetscDraw *sdraw)
694: {
696:   PetscDraw_X    *sXwin = (PetscDraw_X*)(*sdraw)->data;

699:   XFreeGC(sXwin->disp,sXwin->gc.set);
700:   XCloseDisplay(sXwin->disp);
701:   PetscDrawDestroy(&(*sdraw)->popup);
702:   PetscFree((*sdraw)->title);
703:   PetscFree((*sdraw)->display);
704:   PetscFree(sXwin->font);
705:   PetscFree(sXwin);
706:   PetscHeaderDestroy(sdraw);
707:   return(0);
708: }

712: PetscErrorCode PetscDrawXGetDisplaySize_Private(const char name[],int *width,int *height)
713: {
714:   Display *display;

717:   display = XOpenDisplay(name);
718:   if (!display) {
719:     *width  = 0;
720:     *height = 0;
721:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Unable to open display on %s\n.  Make sure your COMPUTE NODES are authorized to connect \n\
722:     to this X server and either your DISPLAY variable\n\
723:     is set or you use the -display name option\n",name);
724:   }
725:   *width  = DisplayWidth(display,0);
726:   *height = DisplayHeight(display,0);
727:   XCloseDisplay(display);
728:   return(0);
729: }

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

734:    Options Database Keys:
735: +  -display <display> - sets the display to use
736: .  -x_virtual - forces use of a X virtual display Xvfb that will not display anything but -draw_save will still work. Xvfb is automatically
737:                 started up in PetscSetDisplay() with this option
738: .  -geometry x,y,w,h - set location and size in pixels
739: -  -draw_size w,h - percentage of screeen, either 1, .5, .3, .25

741:    Level: beginner

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

745: M*/

749: PETSC_EXTERN PetscErrorCode PetscDrawCreate_X(PetscDraw draw)
750: {
751:   PetscDraw_X    *Xwin;
753:   PetscMPIInt    rank;
754:   PetscInt       xywh[4],osize = 4,nsizes=2;
755:   int            x          = draw->x,y = draw->y,w = draw->w,h = draw->h;
756:   static int     xavailable = 0,yavailable = 0,xmax = 0,ymax = 0,ybottom = 0;
757:   PetscBool      flg        = PETSC_FALSE,set;
758:   PetscReal      sizes[2] = {.3,.3};

761:   if (!draw->display) {
762:     PetscMalloc1(256,&draw->display);
763:     PetscGetDisplay(draw->display,256);
764:   }

766:   /*
767:       Initialize the display size
768:   */
769:   if (!xmax) {
770:     PetscDrawXGetDisplaySize_Private(draw->display,&xmax,&ymax);
771:     /* if some processors fail on this and others succed then this is a problem ! */
772:     if (ierr) {
773:       (*PetscErrorPrintf)("PETSc unable to use X windows\nproceeding without graphics\n");
774:       PetscDrawSetType(draw,PETSC_DRAW_NULL);
775:       return(0);
776:     }
777:   }

779:   PetscOptionsGetRealArray(NULL,"-draw_size",sizes,&nsizes,&set);
780:   if (set) {
781:     if (sizes[0] == 1.0)      w = PETSC_DRAW_FULL_SIZE;
782:     else if (sizes[0] == .5)  w = PETSC_DRAW_HALF_SIZE;
783:     else if (sizes[0] == .3)  w = PETSC_DRAW_THIRD_SIZE;
784:     else if (sizes[0] == .25) w = PETSC_DRAW_QUARTER_SIZE;
785:     if (sizes[1] == 1.0)      h = PETSC_DRAW_FULL_SIZE;
786:     else if (sizes[1] == .5)  h = PETSC_DRAW_HALF_SIZE;
787:     else if (sizes[1] == .3)  h = PETSC_DRAW_THIRD_SIZE;
788:     else if (sizes[1] == .25) h = PETSC_DRAW_QUARTER_SIZE;
789:   }

791:   if (w == PETSC_DECIDE) w = draw->w = 300;
792:   if (h == PETSC_DECIDE) h = draw->h = 300;
793:   switch (w) {
794:   case PETSC_DRAW_FULL_SIZE: w = draw->w = xmax - 10; break;
795:   case PETSC_DRAW_HALF_SIZE: w = draw->w = (xmax - 20)/2; break;
796:   case PETSC_DRAW_THIRD_SIZE: w = draw->w = (xmax - 30)/3; break;
797:   case PETSC_DRAW_QUARTER_SIZE: w = draw->w = (xmax - 40)/4; break;
798:   }
799:   switch (h) {
800:   case PETSC_DRAW_FULL_SIZE: h = draw->h = ymax - 10; break;
801:   case PETSC_DRAW_HALF_SIZE: h = draw->h = (ymax - 20)/2; break;
802:   case PETSC_DRAW_THIRD_SIZE: h = draw->h = (ymax - 30)/3; break;
803:   case PETSC_DRAW_QUARTER_SIZE: h = draw->h = (ymax - 40)/4; break;
804:   }

806:   /* allow user to set location and size of window */
807:   xywh[0] = x; xywh[1] = y; xywh[2] = w; xywh[3] = h;
808:   PetscOptionsGetIntArray(NULL,"-geometry",xywh,&osize,NULL);

810:   x = (int) xywh[0]; y = (int) xywh[1]; w = (int) xywh[2]; h = (int) xywh[3];


813:   if (draw->x == PETSC_DECIDE || draw->y == PETSC_DECIDE) {
814:     /*
815:        PETSc tries to place windows starting in the upper left corner and
816:        moving across to the right.

818:               --------------------------------------------
819:               |  Region used so far +xavailable,yavailable |
820:               |                     +                      |
821:               |                     +                      |
822:               |++++++++++++++++++++++ybottom               |
823:               |                                            |
824:               |                                            |
825:               |--------------------------------------------|
826:     */
827:     /*  First: can we add it to the right? */
828:     if (xavailable+w+10 <= xmax) {
829:       x       = xavailable;
830:       y       = yavailable;
831:       ybottom = PetscMax(ybottom,y + h + 30);
832:     } else {
833:       /* No, so add it below on the left */
834:       xavailable = 0;
835:       x          = 0;
836:       yavailable = ybottom;
837:       y          = ybottom;
838:       ybottom    = ybottom + h + 30;
839:     }
840:   }
841:   /* update available region */
842:   xavailable = PetscMax(xavailable,x + w + 10);
843:   if (xavailable >= xmax) {
844:     xavailable = 0;
845:     yavailable = yavailable + h + 30;
846:     ybottom    = yavailable;
847:   }
848:   if (yavailable >= ymax) {
849:     y          = 0;
850:     yavailable = 0;
851:     ybottom    = 0;
852:   }

854:   PetscMemcpy(draw->ops,&DvOps,sizeof(DvOps));

856:   /* actually create and open the window */
857:   PetscNew(&Xwin);
858:   PetscLogObjectMemory((PetscObject)draw,sizeof(PetscDraw_X));
859:   MPI_Comm_rank(PetscObjectComm((PetscObject)draw),&rank);

861:   if (!rank) {
862:     if (x < 0 || y < 0)   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative corner of window");
863:     if (w <= 0 || h <= 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative window width or height");
864:     PetscDrawXiQuickWindow(Xwin,draw->display,draw->title,x,y,w,h);
865:     MPI_Bcast(&Xwin->win,1,MPI_UNSIGNED_LONG,0,PetscObjectComm((PetscObject)draw));
866:   } else {
867:     unsigned long win = 0;
868:     MPI_Bcast(&win,1,MPI_UNSIGNED_LONG,0,PetscObjectComm((PetscObject)draw));
869:     PetscDrawXiQuickWindowFromWindow(Xwin,draw->display,win);
870:   }

872:   Xwin->x    = x;
873:   Xwin->y    = y;
874:   Xwin->w    = w;
875:   Xwin->h    = h;
876:   draw->data = (void*)Xwin;

878:   /*
879:     Need barrier here so processor 0 does not destroy the window before other
880:     processors have completed PetscDrawXiQuickWindow()
881:   */
882:   PetscDrawSynchronizedFlush(draw);

884:   flg  = PETSC_TRUE;
885:   PetscOptionsGetBool(NULL,"-draw_double_buffer",&flg,NULL);
886:   if (flg) {
887:     PetscDrawSetDoubleBuffer(draw);
888:   }
889:   return(0);
890: }