Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TCanvas.cxx
Go to the documentation of this file.
1// @(#)root/gpad:$Id$
2// Author: Rene Brun 12/12/94
3
4/*************************************************************************
5 * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
6 * All rights reserved. *
7 * *
8 * For the licensing terms see $ROOTSYS/LICENSE. *
9 * For the list of contributors see $ROOTSYS/README/CREDITS. *
10 *************************************************************************/
11
12#include <cstring>
13#include <cstdlib>
14#include <iostream>
15#include <fstream>
16
17#include "TROOT.h"
18#include "TBuffer.h"
19#include "TCanvas.h"
20#include "TCanvasImp.h"
21#include "TDatime.h"
22#include "TClass.h"
23#include "TStyle.h"
24#include "TBox.h"
25#include "TCanvasImp.h"
26#include "TDialogCanvas.h"
27#include "TGuiFactory.h"
28#include "TEnv.h"
29#include "TError.h"
30#include "TContextMenu.h"
31#include "TControlBar.h"
32#include "TInterpreter.h"
33#include "TApplication.h"
34#include "TColor.h"
35#include "TObjArray.h"
36#include "TVirtualPadEditor.h"
37#include "TVirtualViewer3D.h"
38#include "TPadPainter.h"
39#include "TVirtualGL.h"
40#include "TVirtualPS.h"
41#include "TVirtualX.h"
42#include "TAxis.h"
43#include "TH1.h"
44#include "TGraph.h"
45#include "TMath.h"
46#include "TView.h"
47#include "strlcpy.h"
48#include "snprintf.h"
49
50#include "TVirtualMutex.h"
51
53public:
56
57//*-*x16 macros/layout_canvas
58
60
62
64
65
66/** \class TCanvas
67\ingroup gpad
68
69The Canvas class.
70
71A Canvas is an area mapped to a window directly under the control of the display
72manager. A ROOT session may have several canvases open at any given time.
73
74A Canvas may be subdivided into independent graphical areas: the __Pads__.
75A canvas has a default pad which has the name of the canvas itself.
76An example of a Canvas layout is sketched in the picture below.
77
78\image html gpad_canvas.png
79
80This canvas contains two pads named P1 and P2. Both Canvas, P1 and P2 can be
81moved, grown, shrunk using the normal rules of the Display manager.
82
83Once objects have been drawn in a canvas, they can be edited/moved by pointing
84directly to them. The cursor shape is changed to suggest the type of action that
85one can do on this object. Clicking with the right mouse button on an object
86pops-up a contextmenu with a complete list of actions possible on this object.
87
88A graphical editor may be started from the canvas "View" menu under the menu
89entry "Toolbar".
90
91An interactive HELP is available by clicking on the HELP button at the top right
92of the canvas. It gives a short explanation about the canvas' menus.
93
94A canvas may be automatically divided into pads via `TPad::Divide`.
95
96At creation time, no matter if in interactive or batch mode, the canvas size
97defines the size of the canvas window (including the size of the window
98manager's decoration). To define precisely the graphics area size of a canvas in
99the interactive mode, the following four lines of code should be used:
100~~~ {.cpp}
101 {
102 Double_t w = 600;
103 Double_t h = 600;
104 auto c = new TCanvas("c", "c", w, h);
105 c->SetWindowSize(w + (w - c->GetWw()), h + (h - c->GetWh()));
106 }
107~~~
108and in the batch mode simply do:
109~~~ {.cpp}
110 c->SetCanvasSize(w,h);
111~~~
112
113If the canvas size this exceed the window size, scroll bars will be added to the canvas
114This allows to display very large canvases (even bigger than the screen size). The
115Following example shows how to proceed.
116~~~ {.cpp}
117 {
118 auto c = new TCanvas("c","c");
119 c->SetCanvasSize(1500, 1500);
120 c->SetWindowSize(500, 500);
121 }
122~~~
123*/
124
125////////////////////////////////////////////////////////////////////////////////
126/// Canvas default constructor.
127
128TCanvas::TCanvas(Bool_t build) : TPad(), fDoubleBuffer(0)
129{
130 fPainter = 0;
131 fWindowTopX = 0;
132 fWindowTopY = 0;
133 fWindowWidth = 0;
134 fWindowHeight = 0;
135 fCw = 0;
136 fCh = 0;
137 fXsizeUser = 0;
138 fYsizeUser = 0;
139 fXsizeReal = kDefaultCanvasSize;
140 fYsizeReal = kDefaultCanvasSize;
141 fHighLightColor = gEnv->GetValue("Canvas.HighLightColor", kRed);
142 fEvent = -1;
143 fEventX = -1;
144 fEventY = -1;
145 fSelectedX = 0;
146 fSelectedY = 0;
147 fRetained = kTRUE;
148 fDrawn = kFALSE;
149 fSelected = 0;
150 fClickSelected = 0;
151 fSelectedPad = 0;
152 fClickSelectedPad = 0;
153 fPadSave = 0;
154 fCanvasImp = 0;
155 fContextMenu = 0;
156
157 fUseGL = gStyle->GetCanvasPreferGL();
158
159 if (!build || TClass::IsCallingNew() != TClass::kRealNew) {
160 Constructor();
161 } else {
162 const char *defcanvas = gROOT->GetDefCanvasName();
163 char *cdef;
164
165 auto lc = (TList*)gROOT->GetListOfCanvases();
166 if (lc->FindObject(defcanvas)) {
167 Int_t n = lc->GetSize()+1;
168 while (lc->FindObject(Form("%s_n%d",defcanvas,n))) n++;
169 cdef = StrDup(Form("%s_n%d",defcanvas,n));
170 } else {
171 cdef = StrDup(Form("%s",defcanvas));
172 }
173 Constructor(cdef, cdef, 1);
174 delete [] cdef;
175 }
176}
177
178////////////////////////////////////////////////////////////////////////////////
179/// Canvas default constructor
180
182{
183 if (gThreadXAR) {
184 void *arr[2];
185 arr[1] = this;
186 if ((*gThreadXAR)("CANV", 2, arr, 0)) return;
187 }
188
189 fCanvas = 0;
190 fCanvasID = -1;
191 fCanvasImp = 0;
192 fBatch = kTRUE;
194
195 fContextMenu = 0;
196 fSelected = 0;
197 fClickSelected = 0;
198 fSelectedPad = 0;
200 fPadSave = 0;
204}
205
206////////////////////////////////////////////////////////////////////////////////
207/// Create an embedded canvas, i.e. a canvas that is in a TGCanvas widget
208/// which is placed in a TGFrame. This ctor is only called via the
209/// TRootEmbeddedCanvas class.
210///
211/// If "name" starts with "gl" the canvas is ready to receive GL output.
212
213TCanvas::TCanvas(const char *name, Int_t ww, Int_t wh, Int_t winid) : TPad(), fDoubleBuffer(0)
214{
215 fCanvasImp = 0;
216 fPainter = 0;
217 Init();
218
219 fCanvasID = winid;
220 fWindowTopX = 0;
221 fWindowTopY = 0;
222 fWindowWidth = ww;
223 fWindowHeight = wh;
224 fCw = ww + 4;
225 fCh = wh +28;
226 fBatch = kFALSE;
228
229 //This is a very special ctor. A window exists already!
230 //Can create painter now.
232
233 if (fUseGL) {
234 fGLDevice = gGLManager->CreateGLContext(winid);
235 if (fGLDevice == -1)
236 fUseGL = kFALSE;
237 }
238
240 if (!fCanvasImp) return;
241
243 SetName(name);
244 Build();
245}
246
247////////////////////////////////////////////////////////////////////////////////
248/// Create a new canvas with a predefined size form.
249/// If form < 0 the menubar is not shown.
250///
251/// - form = 1 700x500 at 10,10 (set by TStyle::SetCanvasDefH,W,X,Y)
252/// - form = 2 500x500 at 20,20
253/// - form = 3 500x500 at 30,30
254/// - form = 4 500x500 at 40,40
255/// - form = 5 500x500 at 50,50
256///
257/// If "name" starts with "gl" the canvas is ready to receive GL output.
258
259TCanvas::TCanvas(const char *name, const char *title, Int_t form) : TPad(), fDoubleBuffer(0)
260{
261 fPainter = 0;
263
264 Constructor(name, title, form);
265}
266
267////////////////////////////////////////////////////////////////////////////////
268/// Create a new canvas with a predefined size form.
269/// If form < 0 the menubar is not shown.
270///
271/// - form = 1 700x500 at 10,10 (set by TStyle::SetCanvasDefH,W,X,Y)
272/// - form = 2 500x500 at 20,20
273/// - form = 3 500x500 at 30,30
274/// - form = 4 500x500 at 40,40
275/// - form = 5 500x500 at 50,50
276
277void TCanvas::Constructor(const char *name, const char *title, Int_t form)
278{
279 if (gThreadXAR) {
280 void *arr[6];
281 static Int_t ww = 500;
282 static Int_t wh = 500;
283 arr[1] = this; arr[2] = (void*)name; arr[3] = (void*)title; arr[4] =&ww; arr[5] = &wh;
284 if ((*gThreadXAR)("CANV", 6, arr, 0)) return;
285 }
286
287 Init();
288 SetBit(kMenuBar,1);
289 if (form < 0) {
290 form = -form;
291 SetBit(kMenuBar,0);
292 }
293
294 fCanvas = this;
295
296 fCanvasID = -1;
297 TCanvas *old = (TCanvas*)gROOT->GetListOfCanvases()->FindObject(name);
298 if (old && old->IsOnHeap()) {
299 Warning("Constructor","Deleting canvas with same name: %s",name);
300 delete old;
301 }
302 if (gROOT->IsBatch()) { //We are in Batch mode
304 if (form == 1) {
307 } else {
308 fWindowWidth = 500;
309 fWindowHeight = 500;
310 }
314 if (!fCanvasImp) return;
315 fBatch = kTRUE;
316 } else { //normal mode with a screen window
318 if (form < 1 || form > 5) form = 1;
319 if (form == 1) {
320 UInt_t uh = UInt_t(cx*gStyle->GetCanvasDefH());
321 UInt_t uw = UInt_t(cx*gStyle->GetCanvasDefW());
322 Int_t ux = Int_t(cx*gStyle->GetCanvasDefX());
323 Int_t uy = Int_t(cx*gStyle->GetCanvasDefY());
324 fCanvasImp = gGuiFactory->CreateCanvasImp(this, name, ux, uy, uw, uh);
325 }
326 fCw = 500;
327 fCh = 500;
328 if (form == 2) fCanvasImp = gGuiFactory->CreateCanvasImp(this, name, 20, 20, UInt_t(cx*500), UInt_t(cx*500));
329 if (form == 3) fCanvasImp = gGuiFactory->CreateCanvasImp(this, name, 30, 30, UInt_t(cx*500), UInt_t(cx*500));
330 if (form == 4) fCanvasImp = gGuiFactory->CreateCanvasImp(this, name, 40, 40, UInt_t(cx*500), UInt_t(cx*500));
331 if (form == 5) fCanvasImp = gGuiFactory->CreateCanvasImp(this, name, 50, 50, UInt_t(cx*500), UInt_t(cx*500));
332 if (!fCanvasImp) return;
333
334 if (!gROOT->IsBatch() && fCanvasID == -1)
336
338 fBatch = kFALSE;
339 }
340
342
343 SetName(name);
344 SetTitle(title); // requires fCanvasImp set
345 Build();
346
347 // Popup canvas
348 fCanvasImp->Show();
349}
350
351////////////////////////////////////////////////////////////////////////////////
352/// Create a new canvas at a random position.
353///
354/// \param[in] name canvas name
355/// \param[in] title canvas title
356/// \param[in] ww is the canvas size in pixels along X
357/// (if ww < 0 the menubar is not shown)
358/// \param[in] wh is the canvas size in pixels along Y
359///
360/// If "name" starts with "gl" the canvas is ready to receive GL output.
361
362TCanvas::TCanvas(const char *name, const char *title, Int_t ww, Int_t wh) : TPad(), fDoubleBuffer(0)
363{
364 fPainter = 0;
366
367 Constructor(name, title, ww, wh);
368}
369
370////////////////////////////////////////////////////////////////////////////////
371/// Create a new canvas at a random position.
372///
373/// \param[in] name canvas name
374/// \param[in] title canvas title
375/// \param[in] ww is the canvas size in pixels along X
376/// (if ww < 0 the menubar is not shown)
377/// \param[in] wh is the canvas size in pixels along Y
378
379void TCanvas::Constructor(const char *name, const char *title, Int_t ww, Int_t wh)
380{
381 if (gThreadXAR) {
382 void *arr[6];
383 arr[1] = this; arr[2] = (void*)name; arr[3] = (void*)title; arr[4] =&ww; arr[5] = &wh;
384 if ((*gThreadXAR)("CANV", 6, arr, 0)) return;
385 }
386
387 Init();
388 SetBit(kMenuBar,1);
389 if (ww < 0) {
390 ww = -ww;
391 SetBit(kMenuBar,0);
392 }
393 if (wh <= 0) {
394 Error("Constructor", "Invalid canvas height: %d",wh);
395 return;
396 }
397 fCw = ww;
398 fCh = wh;
399 fCanvasID = -1;
400 TCanvas *old = (TCanvas*)gROOT->GetListOfCanvases()->FindObject(name);
401 if (old && old->IsOnHeap()) {
402 Warning("Constructor","Deleting canvas with same name: %s",name);
403 delete old;
404 }
405 if (gROOT->IsBatch()) { //We are in Batch mode
407 fWindowWidth = ww;
408 fWindowHeight = wh;
409 fCw = ww;
410 fCh = wh;
412 if (!fCanvasImp) return;
413 fBatch = kTRUE;
414 } else {
416 fCanvasImp = gGuiFactory->CreateCanvasImp(this, name, UInt_t(cx*ww), UInt_t(cx*wh));
417 if (!fCanvasImp) return;
418
419 if (!gROOT->IsBatch() && fCanvasID == -1)
421
423 fBatch = kFALSE;
424 }
425
427
428 SetName(name);
429 SetTitle(title); // requires fCanvasImp set
430 Build();
431
432 // Popup canvas
433 fCanvasImp->Show();
434}
435
436////////////////////////////////////////////////////////////////////////////////
437/// Create a new canvas.
438///
439/// \param[in] name canvas name
440/// \param[in] title canvas title
441/// \param[in] wtopx,wtopy are the pixel coordinates of the top left corner of
442/// the canvas (if wtopx < 0) the menubar is not shown)
443/// \param[in] ww is the canvas size in pixels along X
444/// \param[in] wh is the canvas size in pixels along Y
445///
446/// If "name" starts with "gl" the canvas is ready to receive GL output.
447
448TCanvas::TCanvas(const char *name, const char *title, Int_t wtopx, Int_t wtopy, Int_t ww, Int_t wh)
449 : TPad(), fDoubleBuffer(0)
450{
451 fPainter = 0;
453
454 Constructor(name, title, wtopx, wtopy, ww, wh);
455}
456
457////////////////////////////////////////////////////////////////////////////////
458/// Create a new canvas.
459///
460/// \param[in] name canvas name
461/// \param[in] title canvas title
462/// \param[in] wtopx,wtopy are the pixel coordinates of the top left corner of
463/// the canvas (if wtopx < 0) the menubar is not shown)
464/// \param[in] ww is the canvas size in pixels along X
465/// \param[in] wh is the canvas size in pixels along Y
466
467void TCanvas::Constructor(const char *name, const char *title, Int_t wtopx,
468 Int_t wtopy, Int_t ww, Int_t wh)
469{
470 if (gThreadXAR) {
471 void *arr[8];
472 arr[1] = this; arr[2] = (void*)name; arr[3] = (void*)title;
473 arr[4] = &wtopx; arr[5] = &wtopy; arr[6] = &ww; arr[7] = &wh;
474 if ((*gThreadXAR)("CANV", 8, arr, 0)) return;
475 }
476
477 Init();
478 SetBit(kMenuBar,1);
479 if (wtopx < 0) {
480 wtopx = -wtopx;
481 SetBit(kMenuBar,0);
482 }
483 fCw = ww;
484 fCh = wh;
485 fCanvasID = -1;
486 TCanvas *old = (TCanvas*)gROOT->GetListOfCanvases()->FindObject(name);
487 if (old && old->IsOnHeap()) {
488 Warning("Constructor","Deleting canvas with same name: %s",name);
489 delete old;
490 }
491 if (gROOT->IsBatch()) { //We are in Batch mode
493 fWindowWidth = ww;
494 fWindowHeight = wh;
495 fCw = ww;
496 fCh = wh;
498 if (!fCanvasImp) return;
499 fBatch = kTRUE;
500 } else { //normal mode with a screen window
502 fCanvasImp = gGuiFactory->CreateCanvasImp(this, name, Int_t(cx*wtopx), Int_t(cx*wtopy), UInt_t(cx*ww), UInt_t(cx*wh));
503 if (!fCanvasImp) return;
504
505 if (!gROOT->IsBatch() && fCanvasID == -1)
507
509 fBatch = kFALSE;
510 }
511
513
514 SetName(name);
515 SetTitle(title); // requires fCanvasImp set
516 Build();
517
518 // Popup canvas
519 fCanvasImp->Show();
520}
521
522////////////////////////////////////////////////////////////////////////////////
523/// Initialize the TCanvas members. Called by all constructors.
524
526{
527 // Make sure the application environment exists. It is need for graphics
528 // (colors are initialized in the TApplication ctor).
529 if (!gApplication)
531
532 // Load and initialize graphics libraries if
533 // TApplication::NeedGraphicsLibs() has been called by a
534 // library static initializer.
535 if (gApplication)
537
538 // Get some default from .rootrc. Used in fCanvasImp->InitWindow().
539 fHighLightColor = gEnv->GetValue("Canvas.HighLightColor", kRed);
540 SetBit(kMoveOpaque, gEnv->GetValue("Canvas.MoveOpaque", 0));
541 SetBit(kResizeOpaque, gEnv->GetValue("Canvas.ResizeOpaque", 0));
542 if (gEnv->GetValue("Canvas.ShowEventStatus", kFALSE)) SetBit(kShowEventStatus);
543 if (gEnv->GetValue("Canvas.ShowToolTips", kFALSE)) SetBit(kShowToolTips);
544 if (gEnv->GetValue("Canvas.ShowToolBar", kFALSE)) SetBit(kShowToolBar);
545 if (gEnv->GetValue("Canvas.ShowEditor", kFALSE)) SetBit(kShowEditor);
546 if (gEnv->GetValue("Canvas.AutoExec", kTRUE)) SetBit(kAutoExec);
547
548 // Fill canvas ROOT data structure
549 fXsizeUser = 0;
550 fYsizeUser = 0;
553
554 fDISPLAY = "$DISPLAY";
557 fSelected = 0;
558 fClickSelected = 0;
559 fSelectedX = 0;
560 fSelectedY = 0;
561 fSelectedPad = 0;
563 fPadSave = 0;
564 fEvent = -1;
565 fEventX = -1;
566 fEventY = -1;
567 fContextMenu = 0;
568 fDrawn = kFALSE;
569}
570
571////////////////////////////////////////////////////////////////////////////////
572/// Build a canvas. Called by all constructors.
573
575{
576 // Get window identifier
577 if (fCanvasID == -1 && fCanvasImp)
579 if (fCanvasID == -1) return;
580
581 if (fCw !=0 && fCh !=0) {
584 }
585
586 // Set Pad parameters
587 gPad = this;
588 fCanvas = this;
589 fMother = (TPad*)gPad;
590
591 if (IsBatch()) {
592 // Make sure that batch interactive canvas sizes are the same
593 fCw -= 4;
594 fCh -= 28;
595 } else if (IsWeb()) {
596 // mark canvas as batch - avoid virtualx in many places
598 } else {
599 //normal mode with a screen window
600 // Set default physical canvas attributes
601 //Should be done via gVirtualX, not via fPainter (at least now). No changes here.
602 gVirtualX->SelectWindow(fCanvasID);
603 gVirtualX->SetFillColor(1); //Set color index for fill area
604 gVirtualX->SetLineColor(1); //Set color index for lines
605 gVirtualX->SetMarkerColor(1); //Set color index for markers
606 gVirtualX->SetTextColor(1); //Set color index for text
607 // Clear workstation
608 gVirtualX->ClearWindow();
609
610 // Set Double Buffer on by default
612
613 // Get effective window parameters (with borders and menubar)
616
617 // Get effective canvas parameters without borders
618 Int_t dum1, dum2;
619 gVirtualX->GetGeometry(fCanvasID, dum1, dum2, fCw, fCh);
620
621 fContextMenu = new TContextMenu("ContextMenu");
622 }
623
624 gROOT->GetListOfCanvases()->Add(this);
625
626 if (!fPrimitives) {
627 fPrimitives = new TList;
629 SetFillStyle(1001);
641 fBorderMode=gStyle->GetCanvasBorderMode(); // do not call SetBorderMode (function redefined in TCanvas)
642 SetPad(0, 0, 1, 1);
643 Range(0, 0, 1, 1); //pad range is set by default to [0,1] in x and y
644
646 if (vpp) vpp->SelectDrawable(fPixmapID);//gVirtualX->SelectPixmap(fPixmapID); //pixmap must be selected
647 PaintBorder(GetFillColor(), kTRUE); //paint background
648 }
649
650 // transient canvases have typically no menubar and should not get
651 // by default the event status bar (if set by default)
652 if (TestBit(kMenuBar) && fCanvasImp) {
654 // ... and toolbar + editor
658 }
659}
660
661////////////////////////////////////////////////////////////////////////////////
662/// Canvas destructor
663
665{
666 Destructor();
667}
668
669////////////////////////////////////////////////////////////////////////////////
670/// Browse.
671
673{
674 Draw();
675 cd();
677}
678
679////////////////////////////////////////////////////////////////////////////////
680/// Actual canvas destructor.
681
683{
684 if (gThreadXAR) {
685 void *arr[2];
686 arr[1] = this;
687 if ((*gThreadXAR)("CDEL", 2, arr, 0)) return;
688 }
689
690 if (ROOT::Detail::HasBeenDeleted(this)) return;
691
693 if (!gPad) return;
694
695 Close();
696
697 //If not yet (batch mode?).
699}
700
701////////////////////////////////////////////////////////////////////////////////
702/// Set current canvas & pad. Returns the new current pad,
703/// or 0 in case of failure.
704/// See TPad::cd() for an explanation of the parameter.
705
707{
708 if (fCanvasID == -1) return 0;
709
710 TPad::cd(subpadnumber);
711
712 // in case doublebuffer is off, draw directly onto display window
713 if (!IsBatch()) {
714 if (!fDoubleBuffer)
715 gVirtualX->SelectWindow(fCanvasID);//Ok, does not matter for glpad.
716 }
717 return gPad;
718}
719
720////////////////////////////////////////////////////////////////////////////////
721/// Remove all primitives from the canvas.
722/// If option "D" is specified, direct sub-pads are cleared but not deleted.
723/// This option is not recursive, i.e. pads in direct sub-pads are deleted.
724
726{
727 if (fCanvasID == -1) return;
728
730
731 TString opt = option;
732 opt.ToLower();
733 if (opt.Contains("d")) {
734 // clear subpads, but do not delete pads in case the canvas
735 // has been divided (note: option "D" is propagated so could cause
736 // conflicts for primitives using option "D" for something else)
737 if (fPrimitives) {
738 TIter next(fPrimitives);
739 TObject *obj;
740 while ((obj=next())) {
741 obj->Clear(option);
742 }
743 }
744 } else {
745 //default, clear everything in the canvas. Subpads are deleted
746 TPad::Clear(option); //Remove primitives from pad
747 }
748
749 fSelected = 0;
750 fClickSelected = 0;
751 fSelectedPad = 0;
753}
754
755////////////////////////////////////////////////////////////////////////////////
756/// Emit pad Cleared signal.
757
759{
760 Emit("Cleared(TVirtualPad*)", (Longptr_t)pad);
761}
762
763////////////////////////////////////////////////////////////////////////////////
764/// Emit Closed signal.
765
767{
768 Emit("Closed()");
769}
770
771////////////////////////////////////////////////////////////////////////////////
772/// Close canvas.
773///
774/// Delete window/pads data structure
775
777{
778 TPad *padsave = (TPad*)gPad;
779 TCanvas *cansave = 0;
780 if (padsave) cansave = (TCanvas*)gPad->GetCanvas();
781
782 if (fCanvasID != -1) {
783
784 if ((!gROOT->IsLineProcessing()) && (!gVirtualX->IsCmdThread())) {
785 gInterpreter->Execute(this, IsA(), "Close", option);
786 return;
787 }
788
790
792
793 cd();
794 TPad::Close(option);
795
796 if (!IsBatch()) {
797 gVirtualX->SelectWindow(fCanvasID); //select current canvas
798
800
801 if (fCanvasImp)
802 fCanvasImp->Close();
803 }
804 fCanvasID = -1;
805 fBatch = kTRUE;
806
807 gROOT->GetListOfCanvases()->Remove(this);
808
809 // Close actual window on screen
811 }
812
813 if (cansave == this) {
814 gPad = (TCanvas *) gROOT->GetListOfCanvases()->First();
815 } else {
816 gPad = padsave;
817 }
818
819 Closed();
820}
821
822////////////////////////////////////////////////////////////////////////////////
823/// Copy the canvas pixmap of the pad to the canvas.
824
826{
827 if (!IsBatch()) {
828 CopyPixmap();
830 }
831}
832
833////////////////////////////////////////////////////////////////////////////////
834/// Draw a canvas.
835/// If a canvas with the name is already on the screen, the canvas is repainted.
836/// This function is useful when a canvas object has been saved in a Root file.
837/// One can then do:
838/// ~~~ {.cpp}
839/// Root > Tfile f("file.root");
840/// Root > canvas.Draw();
841/// ~~~
842
844{
845 // Load and initialize graphics libraries if
846 // TApplication::NeedGraphicsLibs() has been called by a
847 // library static initializer.
848 if (gApplication)
850
851 fDrawn = kTRUE;
852
853 TCanvas *old = (TCanvas*)gROOT->GetListOfCanvases()->FindObject(GetName());
854 if (old == this) {
855 Paint();
856 return;
857 }
858 if (old) { gROOT->GetListOfCanvases()->Remove(old); delete old;}
859
860 if (fWindowWidth == 0) {
861 if (fCw !=0) fWindowWidth = fCw+4;
862 else fWindowWidth = 800;
863 }
864 if (fWindowHeight == 0) {
865 if (fCh !=0) fWindowHeight = fCh+28;
866 else fWindowHeight = 600;
867 }
868 if (gROOT->IsBatch()) { //We are in Batch mode
870 if (!fCanvasImp) return;
871 fBatch = kTRUE;
872
873 } else { //normal mode with a screen window
876 if (!fCanvasImp) return;
878 }
879 Build();
880 ResizePad();
882 fCanvasImp->Show();
883 Modified();
884}
885
886////////////////////////////////////////////////////////////////////////////////
887/// Draw a clone of this canvas
888/// A new canvas is created that is a clone of this canvas
889
891{
892 TCanvas *newCanvas = (TCanvas*)Clone();
893 newCanvas->SetName();
894
895 newCanvas->Draw(option);
896 newCanvas->Update();
897 return newCanvas;
898}
899
900////////////////////////////////////////////////////////////////////////////////
901/// Draw a clone of this canvas into the current pad
902/// In an interactive session, select the destination/current pad
903/// with the middle mouse button, then point to the canvas area to select
904/// the canvas context menu item DrawClonePad.
905/// Note that the original canvas may have subpads.
906
908{
909 TPad *padsav = (TPad*)gPad;
910 TPad *selpad = (TPad*)gROOT->GetSelectedPad();
911 TPad *pad = padsav;
912 if (pad == this) pad = selpad;
913 if (padsav == 0 || pad == 0 || pad == this) {
914 TCanvas *newCanvas = (TCanvas*)DrawClone();
916 return newCanvas;
917 }
918 if (fCanvasID == -1) {
921 if (!fCanvasImp) return 0;
924 }
925 this->cd();
926 TObject *obj, *clone;
927 //copy pad attributes
928 pad->Range(fX1,fY1,fX2,fY2);
929 pad->SetTickx(GetTickx());
930 pad->SetTicky(GetTicky());
931 pad->SetGridx(GetGridx());
932 pad->SetGridy(GetGridy());
933 pad->SetLogx(GetLogx());
934 pad->SetLogy(GetLogy());
935 pad->SetLogz(GetLogz());
938 TAttLine::Copy((TAttLine&)*pad);
939 TAttFill::Copy((TAttFill&)*pad);
940 TAttPad::Copy((TAttPad&)*pad);
941
942 //copy primitives
944 while ((obj=next())) {
945 pad->cd();
946 clone = obj->Clone();
947 pad->GetListOfPrimitives()->Add(clone,next.GetOption());
948 }
949 pad->ResizePad();
950 pad->Modified();
951 pad->Update();
952 if (padsav) padsav->cd();
953 return 0;
954}
955
956////////////////////////////////////////////////////////////////////////////////
957/// Report name and title of primitive below the cursor.
958///
959/// This function is called when the option "Event Status"
960/// in the canvas menu "Options" is selected.
961
963{
964 const Int_t kTMAX=256;
965 static char atext[kTMAX];
966
967 if (!TestBit(kShowEventStatus) || !selected) return;
968
969 if (!fCanvasImp) return; //this may happen when closing a TAttCanvas
970
971 TVirtualPad* savepad;
972 savepad = gPad;
974
975 fCanvasImp->SetStatusText(selected->GetTitle(),0);
976 fCanvasImp->SetStatusText(selected->GetName(),1);
977 if (event == kKeyPress)
978 snprintf(atext, kTMAX, "%c", (char) px);
979 else
980 snprintf(atext, kTMAX, "%d,%d", px, py);
981 fCanvasImp->SetStatusText(atext,2);
982
983 // Show date/time if TimeDisplay is selected
984 TAxis *xaxis = NULL;
985 if ( selected->InheritsFrom("TH1") )
986 xaxis = ((TH1*)selected)->GetXaxis();
987 else if ( selected->InheritsFrom("TGraph") )
988 xaxis = ((TGraph*)selected)->GetXaxis();
989 else if ( selected->InheritsFrom("TAxis") )
990 xaxis = (TAxis*)selected;
991 if ( xaxis != NULL && xaxis->GetTimeDisplay()) {
992 TString objinfo = selected->GetObjectInfo(px,py);
993 // check if user has overwritten GetObjectInfo and altered
994 // the default text from TObject::GetObjectInfo "x=.. y=.."
995 if (objinfo.Contains("x=") && objinfo.Contains("y=") ) {
996 UInt_t toff = 0;
997 TString time_format(xaxis->GetTimeFormat());
998 // TimeFormat may contain offset: %F2000-01-01 00:00:00
999 Int_t idF = time_format.Index("%F");
1000 if (idF>=0) {
1001 Int_t lnF = time_format.Length();
1002 // minimal check for correct format
1003 if (lnF - idF == 21) {
1004 time_format = time_format(idF+2, lnF);
1005 TDatime dtoff(time_format);
1006 toff = dtoff.Convert();
1007 }
1008 } else {
1009 toff = (UInt_t)gStyle->GetTimeOffset();
1010 }
1011 TDatime dt((UInt_t)gPad->AbsPixeltoX(px) + toff);
1012 snprintf(atext, kTMAX, "%s, y=%g",
1013 dt.AsSQLString(),gPad->AbsPixeltoY(py));
1014 fCanvasImp->SetStatusText(atext,3);
1015 gPad = savepad;
1016 return;
1017 }
1018 }
1019 // default
1020 fCanvasImp->SetStatusText(selected->GetObjectInfo(px,py),3);
1021
1022 gPad = savepad;
1023}
1024
1025////////////////////////////////////////////////////////////////////////////////
1026/// Get editor bar.
1027
1029{
1031}
1032
1033////////////////////////////////////////////////////////////////////////////////
1034/// Embedded a canvas into a TRootEmbeddedCanvas. This method is only called
1035/// via TRootEmbeddedCanvas::AdoptCanvas.
1036
1038{
1039 // If fCanvasImp already exists, no need to go further.
1040 if(fCanvasImp) return;
1041
1042 fCanvasID = winid;
1043 fWindowTopX = 0;
1044 fWindowTopY = 0;
1045 fWindowWidth = ww;
1046 fWindowHeight = wh;
1047 fCw = ww;
1048 fCh = wh;
1049 fBatch = kFALSE;
1050 fUpdating = kFALSE;
1051
1053 if (!fCanvasImp) return;
1054 Build();
1055 Resize();
1056}
1057
1058////////////////////////////////////////////////////////////////////////////////
1059/// Generate kMouseEnter and kMouseLeave events depending on the previously
1060/// selected object and the currently selected object. Does nothing if the
1061/// selected object does not change.
1062
1063void TCanvas::EnterLeave(TPad *prevSelPad, TObject *prevSelObj)
1064{
1065 if (prevSelObj == fSelected) return;
1066
1067 TPad *padsav = (TPad *)gPad;
1068 Int_t sevent = fEvent;
1069
1070 if (prevSelObj) {
1071 gPad = prevSelPad;
1072 prevSelObj->ExecuteEvent(kMouseLeave, fEventX, fEventY);
1074 RunAutoExec();
1075 ProcessedEvent(kMouseLeave, fEventX, fEventY, prevSelObj); // emit signal
1076 }
1077
1079
1080 if (fSelected) {
1083 RunAutoExec();
1085 }
1086
1087 fEvent = sevent;
1088 gPad = padsav;
1089}
1090
1091////////////////////////////////////////////////////////////////////////////////
1092/// Execute action corresponding to one event.
1093///
1094/// This member function must be implemented to realize the action
1095/// corresponding to the mouse click on the object in the canvas
1096///
1097/// Only handle mouse motion events in TCanvas, all other events are
1098/// ignored for the time being
1099
1101{
1102 if (gROOT->GetEditorMode()) {
1104 return;
1105 }
1106
1107 switch (event) {
1108
1109 case kMouseMotion:
1111 break;
1112 }
1113}
1114
1115////////////////////////////////////////////////////////////////////////////////
1116/// Turn rubberband feedback mode on or off.
1117
1119{
1120 if (set) {
1121 SetDoubleBuffer(0); // turn off double buffer mode
1122 gVirtualX->SetDrawMode(TVirtualX::kInvert); // set the drawing mode to XOR mode
1123 } else {
1124 SetDoubleBuffer(1); // turn on double buffer mode
1125 gVirtualX->SetDrawMode(TVirtualX::kCopy); // set drawing mode back to normal (copy) mode
1126 }
1127}
1128
1129////////////////////////////////////////////////////////////////////////////////
1130/// Flush canvas buffers.
1131
1133{
1134 if ((fCanvasID == -1) || IsWeb()) return;
1135
1136 TPad *padsav = (TPad*)gPad;
1137 cd();
1138 if (!IsBatch()) {
1139 if (!UseGL() || fGLDevice == -1) {
1140 gVirtualX->SelectWindow(fCanvasID);
1141 gPad = padsav; //don't do cd() because than also the pixmap is changed
1142 CopyPixmaps();
1143 gVirtualX->UpdateWindow(1);
1144 } else {
1145 TVirtualPS *tvps = gVirtualPS;
1146 gVirtualPS = 0;
1147 gGLManager->MakeCurrent(fGLDevice);
1149 Paint();
1150 if (padsav && padsav->GetCanvas() == this) {
1151 padsav->cd();
1152 padsav->HighLight(padsav->GetHighLightColor());
1153 //cd();
1154 }
1156 gGLManager->Flush(fGLDevice);
1157 gVirtualPS = tvps;
1158 }
1159 }
1160 if (padsav) padsav->cd();
1161}
1162
1163////////////////////////////////////////////////////////////////////////////////
1164/// Force canvas update
1165
1167{
1169}
1170
1171////////////////////////////////////////////////////////////////////////////////
1172/// Force a copy of current style for all objects in canvas.
1173
1175{
1176 if ((!gROOT->IsLineProcessing()) && (!gVirtualX->IsCmdThread())) {
1177 gInterpreter->Execute(this, IsA(), "UseCurrentStyle", "");
1178 return;
1179 }
1180
1182
1184
1185 if (gStyle->IsReading()) {
1189 } else {
1193 }
1194}
1195
1196////////////////////////////////////////////////////////////////////////////////
1197/// Returns current top x position of window on screen.
1198
1200{
1203
1204 return fWindowTopX;
1205}
1206
1207////////////////////////////////////////////////////////////////////////////////
1208/// Returns current top y position of window on screen.
1209
1211{
1214
1215 return fWindowTopY;
1216}
1217
1218////////////////////////////////////////////////////////////////////////////////
1219/// Handle Input Events.
1220///
1221/// Handle input events, like button up/down in current canvas.
1222
1224{
1225 TPad *pad;
1226 TPad *prevSelPad = (TPad*) fSelectedPad;
1227 TObject *prevSelObj = fSelected;
1228
1229 fPadSave = (TPad*)gPad;
1230 cd(); // make sure this canvas is the current canvas
1231
1232 fEvent = event;
1233 fEventX = px;
1234 fEventY = py;
1235
1236 switch (event) {
1237
1238 case kMouseMotion:
1239 // highlight object tracked over
1240 pad = Pick(px, py, prevSelObj);
1241 if (!pad) return;
1242
1243 EnterLeave(prevSelPad, prevSelObj);
1244
1245 gPad = pad; // don't use cd() we will use the current
1246 // canvas via the GetCanvas member and not via
1247 // gPad->GetCanvas
1248
1249 if (fSelected) {
1250 fSelected->ExecuteEvent(event, px, py);
1251 RunAutoExec();
1252 }
1253
1254 break;
1255
1256 case kMouseEnter:
1257 // mouse enters canvas
1259 break;
1260
1261 case kMouseLeave:
1262 // mouse leaves canvas
1263 {
1264 // force popdown of tooltips
1265 TObject *sobj = fSelected;
1266 TPad *spad = fSelectedPad;
1267 fSelected = 0;
1268 fSelectedPad = 0;
1269 EnterLeave(prevSelPad, prevSelObj);
1270 fSelected = sobj;
1271 fSelectedPad = spad;
1273 }
1274 break;
1275
1276 case kButton1Double:
1277 // triggered on the second button down within 350ms and within
1278 // 3x3 pixels of the first button down, button up finishes action
1279
1280 case kButton1Down:
1281 // find pad in which input occurred
1282 pad = Pick(px, py, prevSelObj);
1283 if (!pad) return;
1284
1285 gPad = pad; // don't use cd() because we won't draw in pad
1286 // we will only use its coordinate system
1287
1288 if (fSelected) {
1289 FeedbackMode(kTRUE); // to draw in rubberband mode
1290 fSelected->ExecuteEvent(event, px, py);
1291
1292 RunAutoExec();
1293 }
1294
1295 break;
1296
1297 case kArrowKeyPress:
1298 case kArrowKeyRelease:
1299 case kButton1Motion:
1300 case kButton1ShiftMotion: //8 == kButton1Motion + shift modifier
1301 if (fSelected) {
1303
1304 fSelected->ExecuteEvent(event, px, py);
1305 gVirtualX->Update();
1306 if (fSelected && !fSelected->InheritsFrom(TAxis::Class())) {
1307 Bool_t resize = kFALSE;
1308 if (fSelected->InheritsFrom(TBox::Class()))
1309 resize = ((TBox*)fSelected)->IsBeingResized();
1310 if (fSelected->InheritsFrom(TVirtualPad::Class()))
1311 resize = ((TVirtualPad*)fSelected)->IsBeingResized();
1312
1313 if ((!resize && TestBit(kMoveOpaque)) || (resize && TestBit(kResizeOpaque))) {
1314 gPad = fPadSave;
1315 Update();
1317 }
1318 }
1319
1320 RunAutoExec();
1321 }
1322
1323 break;
1324
1325 case kButton1Up:
1326
1327 if (fSelected) {
1329
1330 fSelected->ExecuteEvent(event, px, py);
1331
1332 RunAutoExec();
1333
1334 if (fPadSave)
1335 gPad = fPadSave;
1336 else {
1337 gPad = this;
1338 fPadSave = this;
1339 }
1340
1341 Update(); // before calling update make sure gPad is reset
1342 }
1343 break;
1344
1345//*-*----------------------------------------------------------------------
1346
1347 case kButton2Down:
1348 // find pad in which input occurred
1349 pad = Pick(px, py, prevSelObj);
1350 if (!pad) return;
1351
1352 gPad = pad; // don't use cd() because we won't draw in pad
1353 // we will only use its coordinate system
1354
1356
1357 if (fSelected) fSelected->Pop(); // pop object to foreground
1358 pad->cd(); // and make its pad the current pad
1359 if (gDebug)
1360 printf("Current Pad: %s / %s\n", pad->GetName(), pad->GetTitle());
1361
1362 // loop over all canvases to make sure that only one pad is highlighted
1363 {
1364 TIter next(gROOT->GetListOfCanvases());
1365 TCanvas *tc;
1366 while ((tc = (TCanvas *)next()))
1367 tc->Update();
1368 }
1369
1370 //if (pad->GetGLDevice() != -1 && fSelected)
1371 // fSelected->ExecuteEvent(event, px, py);
1372
1373 break; // don't want fPadSave->cd() to be executed at the end
1374
1375 case kButton2Motion:
1376 //was empty!
1377 case kButton2Up:
1378 if (fSelected) {
1380
1381 fSelected->ExecuteEvent(event, px, py);
1382 RunAutoExec();
1383 }
1384 break;
1385
1386 case kButton2Double:
1387 break;
1388
1389//*-*----------------------------------------------------------------------
1390
1391 case kButton3Down:
1392 // popup context menu
1393 pad = Pick(px, py, prevSelObj);
1394 if (!pad) return;
1395
1397
1400 fContextMenu->Popup(px, py, fSelected, this, pad);
1401
1402 break;
1403
1404 case kButton3Motion:
1405 break;
1406
1407 case kButton3Up:
1409 break;
1410
1411 case kButton3Double:
1412 break;
1413
1414 case kKeyPress:
1415 if (!fSelectedPad || !fSelected) return;
1416 gPad = fSelectedPad; // don't use cd() because we won't draw in pad
1417 // we will only use its coordinate system
1418 fSelected->ExecuteEvent(event, px, py);
1419
1420 RunAutoExec();
1421
1422 break;
1423
1424 case kButton1Shift:
1425 // Try to select
1426 pad = Pick(px, py, prevSelObj);
1427
1428 if (!pad) return;
1429
1430 EnterLeave(prevSelPad, prevSelObj);
1431
1432 gPad = pad; // don't use cd() we will use the current
1433 // canvas via the GetCanvas member and not via
1434 // gPad->GetCanvas
1435 if (fSelected) {
1436 fSelected->ExecuteEvent(event, px, py);
1437 RunAutoExec();
1438 }
1439 break;
1440
1441 case kWheelUp:
1442 case kWheelDown:
1443 pad = Pick(px, py, prevSelObj);
1444 if (!pad) return;
1445
1446 gPad = pad;
1447 if (fSelected)
1448 fSelected->ExecuteEvent(event, px, py);
1449 break;
1450
1451 default:
1452 break;
1453 }
1454
1455 if (fPadSave && event != kButton2Down)
1456 fPadSave->cd();
1457
1458 if (event != kMouseLeave) { // signal was already emitted for this event
1459 ProcessedEvent(event, px, py, fSelected); // emit signal
1461 }
1462}
1463
1464////////////////////////////////////////////////////////////////////////////////
1465/// Iconify canvas
1466
1468{
1469 if (fCanvasImp)
1471}
1472
1473////////////////////////////////////////////////////////////////////////////////
1474/// Is folder ?
1475
1477{
1478 return fgIsFolder;
1479}
1480
1481////////////////////////////////////////////////////////////////////////////////
1482/// Is web canvas
1483
1485{
1486 return fCanvasImp ? fCanvasImp->IsWeb() : kFALSE;
1487}
1488
1489////////////////////////////////////////////////////////////////////////////////
1490/// List all pads.
1491
1492void TCanvas::ls(Option_t *option) const
1493{
1495 std::cout <<"Canvas Name=" <<GetName()<<" Title="<<GetTitle()<<" Option="<<option<<std::endl;
1497 TPad::ls(option);
1499}
1500
1501////////////////////////////////////////////////////////////////////////////////
1502/// Static function to build a default canvas.
1503
1505{
1506 const char *defcanvas = gROOT->GetDefCanvasName();
1507 char *cdef;
1508
1509 auto lc = (TList*)gROOT->GetListOfCanvases();
1510 if (lc->FindObject(defcanvas)) {
1511 Int_t n = lc->GetSize() + 1;
1512 cdef = new char[strlen(defcanvas)+15];
1513 do {
1514 strlcpy(cdef,Form("%s_n%d", defcanvas, n++),strlen(defcanvas)+15);
1515 } while (lc->FindObject(cdef));
1516 } else
1517 cdef = StrDup(Form("%s",defcanvas));
1518
1519 TCanvas *c = new TCanvas(cdef, cdef, 1);
1520
1521 ::Info("TCanvas::MakeDefCanvas"," created default TCanvas with name %s",cdef);
1522 delete [] cdef;
1523 return c;
1524}
1525
1526////////////////////////////////////////////////////////////////////////////////
1527/// Set option to move objects/pads in a canvas.
1528///
1529/// - set = 1 (default) graphics objects are moved in opaque mode
1530/// - set = 0 only the outline of objects is drawn when moving them
1531///
1532/// The option opaque produces the best effect. It requires however a
1533/// a reasonably fast workstation or response time.
1534
1536{
1537 SetBit(kMoveOpaque,set);
1538}
1539
1540////////////////////////////////////////////////////////////////////////////////
1541/// Paint canvas.
1542
1544{
1545 if (fCanvas)
1546 TPad::Paint(option);
1547}
1548
1549////////////////////////////////////////////////////////////////////////////////
1550/// Prepare for pick, call TPad::Pick() and when selected object
1551/// is different from previous then emit Picked() signal.
1552
1553TPad *TCanvas::Pick(Int_t px, Int_t py, TObject *prevSelObj)
1554{
1555 TObjLink *pickobj = 0;
1556
1557 fSelected = 0;
1558 fSelectedOpt = "";
1559 fSelectedPad = 0;
1560
1561 TPad *pad = Pick(px, py, pickobj);
1562 if (!pad) return 0;
1563
1564 if (!pickobj) {
1565 fSelected = pad;
1566 fSelectedOpt = "";
1567 } else {
1568 if (!fSelected) { // can be set via TCanvas::SetSelected()
1569 fSelected = pickobj->GetObject();
1570 fSelectedOpt = pickobj->GetOption();
1571 }
1572 }
1573 fSelectedPad = pad;
1574
1575 if (fSelected != prevSelObj)
1576 Picked(fSelectedPad, fSelected, fEvent); // emit signal
1577
1578 if ((fEvent == kButton1Down) || (fEvent == kButton2Down) || (fEvent == kButton3Down)) {
1579 if (fSelected && !fSelected->InheritsFrom(TView::Class())) {
1582 Selected(fSelectedPad, fSelected, fEvent); // emit signal
1583 fSelectedX = px;
1584 fSelectedY = py;
1585 }
1586 }
1587 return pad;
1588}
1589
1590////////////////////////////////////////////////////////////////////////////////
1591/// Emit Picked() signal.
1592
1594{
1595 Longptr_t args[3];
1596
1597 args[0] = (Longptr_t) pad;
1598 args[1] = (Longptr_t) obj;
1599 args[2] = event;
1600
1601 Emit("Picked(TPad*,TObject*,Int_t)", args);
1602}
1603
1604////////////////////////////////////////////////////////////////////////////////
1605/// Emit Highlighted() signal.
1606///
1607/// - pad is pointer to pad with highlighted histogram or graph
1608/// - obj is pointer to highlighted histogram or graph
1609/// - x is highlighted x bin for 1D histogram or highlighted x-th point for graph
1610/// - y is highlighted y bin for 2D histogram (for 1D histogram or graph not in use)
1611
1613{
1614 Longptr_t args[4];
1615
1616 args[0] = (Longptr_t) pad;
1617 args[1] = (Longptr_t) obj;
1618 args[2] = x;
1619 args[3] = y;
1620
1621 Emit("Highlighted(TVirtualPad*,TObject*,Int_t,Int_t)", args);
1622}
1623
1624////////////////////////////////////////////////////////////////////////////////
1625/// This is "simplification" for function TCanvas::Connect with Highlighted
1626/// signal for specific slot.
1627///
1628/// Slot has to be defined "UserFunction(TVirtualPad *pad, TObject *obj, Int_t x, Int_t y)"
1629/// all parameters of UserFunction are taken from TCanvas::Highlighted
1630
1631void TCanvas::HighlightConnect(const char *slot)
1632{
1633 Connect("Highlighted(TVirtualPad*,TObject*,Int_t,Int_t)", 0, 0, slot);
1634}
1635
1636////////////////////////////////////////////////////////////////////////////////
1637/// Emit Selected() signal.
1638
1640{
1641 Longptr_t args[3];
1642
1643 args[0] = (Longptr_t) pad;
1644 args[1] = (Longptr_t) obj;
1645 args[2] = event;
1646
1647 Emit("Selected(TVirtualPad*,TObject*,Int_t)", args);
1648}
1649
1650////////////////////////////////////////////////////////////////////////////////
1651/// Emit ProcessedEvent() signal.
1652
1654{
1655 Longptr_t args[4];
1656
1657 args[0] = event;
1658 args[1] = x;
1659 args[2] = y;
1660 args[3] = (Longptr_t) obj;
1661
1662 Emit("ProcessedEvent(Int_t,Int_t,Int_t,TObject*)", args);
1663}
1664
1665////////////////////////////////////////////////////////////////////////////////
1666/// Recompute canvas parameters following a X11 Resize.
1667
1669{
1670 if (fCanvasID == -1) return;
1671
1672 if ((!gROOT->IsLineProcessing()) && (!gVirtualX->IsCmdThread())) {
1673 gInterpreter->Execute(this, IsA(), "Resize", "");
1674 return;
1675 }
1676
1678
1679 TPad *padsav = (TPad*)gPad;
1680 cd();
1681
1682 if (!IsBatch()) {
1683 gVirtualX->SelectWindow(fCanvasID); //select current canvas
1684 gVirtualX->ResizeWindow(fCanvasID); //resize canvas and off-screen buffer
1685
1686 // Get effective window parameters including menubar and borders
1689
1690 // Get effective canvas parameters without borders
1691 Int_t dum1, dum2;
1692 gVirtualX->GetGeometry(fCanvasID, dum1, dum2, fCw, fCh);
1693 }
1694
1695 if (fXsizeUser && fYsizeUser) {
1696 UInt_t nwh = fCh;
1697 UInt_t nww = fCw;
1699 if (rxy < 1) {
1700 UInt_t twh = UInt_t(Double_t(fCw)/rxy);
1701 if (twh > fCh)
1702 nww = UInt_t(Double_t(fCh)*rxy);
1703 else
1704 nwh = twh;
1705 if (nww > fCw) {
1706 nww = fCw; nwh = twh;
1707 }
1708 if (nwh > fCh) {
1709 nwh = fCh; nww = UInt_t(Double_t(fCh)/rxy);
1710 }
1711 } else {
1712 UInt_t twh = UInt_t(Double_t(fCw)*rxy);
1713 if (twh > fCh)
1714 nwh = UInt_t(Double_t(fCw)/rxy);
1715 else
1716 nww = twh;
1717 if (nww > fCw) {
1718 nww = fCw; nwh = twh;
1719 }
1720 if (nwh > fCh) {
1721 nwh = fCh; nww = UInt_t(Double_t(fCh)*rxy);
1722 }
1723 }
1724 fCw = nww;
1725 fCh = nwh;
1726 }
1727
1728 if (fCw < fCh) {
1731 }
1732 else {
1735 }
1736
1737//*-*- Loop on all pads to recompute conversion coefficients
1739
1740 if (padsav) padsav->cd();
1741}
1742
1743
1744////////////////////////////////////////////////////////////////////////////////
1745/// Raise canvas window
1746
1748{
1749 if (fCanvasImp)
1751}
1752
1753////////////////////////////////////////////////////////////////////////////////
1754/// Set option to resize objects/pads in a canvas.
1755///
1756/// - set = 1 (default) graphics objects are resized in opaque mode
1757/// - set = 0 only the outline of objects is drawn when resizing them
1758///
1759/// The option opaque produces the best effect. It requires however a
1760/// a reasonably fast workstation or response time.
1761
1763{
1764 SetBit(kResizeOpaque,set);
1765}
1766
1767////////////////////////////////////////////////////////////////////////////////
1768/// Execute the list of TExecs in the current pad.
1769
1771{
1772 if (!TestBit(kAutoExec)) return;
1773 if (!gPad) return;
1774 ((TPad*)gPad)->AutoExec();
1775}
1776
1777
1778////////////////////////////////////////////////////////////////////////////////
1779/// Save primitives in this canvas in C++ macro file with GUI.
1780
1781void TCanvas::SavePrimitive(std::ostream &out, Option_t *option /*= ""*/)
1782{
1783 // Write canvas options (in $TROOT or $TStyle)
1784 if (gStyle->GetOptFit()) {
1785 out<<" gStyle->SetOptFit(1);"<<std::endl;
1786 }
1787 if (!gStyle->GetOptStat()) {
1788 out<<" gStyle->SetOptStat(0);"<<std::endl;
1789 }
1790 if (!gStyle->GetOptTitle()) {
1791 out<<" gStyle->SetOptTitle(0);"<<std::endl;
1792 }
1793 if (gROOT->GetEditHistograms()) {
1794 out<<" gROOT->SetEditHistograms();"<<std::endl;
1795 }
1796 if (GetShowEventStatus()) {
1797 out<<" "<<GetName()<<"->ToggleEventStatus();"<<std::endl;
1798 }
1799 if (GetShowToolTips()) {
1800 out<<" "<<GetName()<<"->ToggleToolTips();"<<std::endl;
1801 }
1802 if (GetShowToolBar()) {
1803 out<<" "<<GetName()<<"->ToggleToolBar();"<<std::endl;
1804 }
1805 if (GetHighLightColor() != 5) {
1806 if (GetHighLightColor() > 228) {
1808 out<<" "<<GetName()<<"->SetHighLightColor(ci);" << std::endl;
1809 } else
1810 out<<" "<<GetName()<<"->SetHighLightColor("<<GetHighLightColor()<<");"<<std::endl;
1811 }
1812
1813 // Now recursively scan all pads of this canvas
1814 cd();
1815 TPad::SavePrimitive(out,option);
1816}
1817
1818////////////////////////////////////////////////////////////////////////////////
1819/// Save primitives in this canvas as a C++ macro file.
1820/// This function loops on all the canvas primitives and for each primitive
1821/// calls the object SavePrimitive function.
1822/// When outputting floating point numbers, the default precision is 7 digits.
1823/// The precision can be changed (via system.rootrc) by changing the value
1824/// of the environment variable "Canvas.SavePrecision"
1825
1826void TCanvas::SaveSource(const char *filename, Option_t *option)
1827{
1828 // reset bit TClass::kClassSaved for all classes
1829 TIter next(gROOT->GetListOfClasses());
1830 TClass *cl;
1831 while((cl = (TClass*)next())) {
1833 }
1834
1835 char quote = '"';
1836 std::ofstream out;
1837 Int_t lenfile = strlen(filename);
1838 char * fname;
1839 char lcname[10];
1840 const char *cname = GetName();
1841 Bool_t invalid = kFALSE;
1842 // if filename is given, open this file, otherwise create a file
1843 // with a name equal to the canvasname.C
1844 if (lenfile) {
1845 fname = (char*)filename;
1846 out.open(fname, std::ios::out);
1847 } else {
1848 Int_t nch = strlen(cname);
1849 if (nch < 10) {
1850 strlcpy(lcname,cname,10);
1851 for (Int_t k=1;k<=nch;k++) {if (lcname[nch-k] == ' ') lcname[nch-k] = 0;}
1852 if (lcname[0] == 0) {invalid = kTRUE; strlcpy(lcname,"c1",10); nch = 2;}
1853 cname = lcname;
1854 }
1855 fname = new char[nch+3];
1856 strlcpy(fname,cname,nch+3);
1857 strncat(fname,".C",3);
1858 out.open(fname, std::ios::out);
1859 }
1860 if (!out.good ()) {
1861 Error("SaveSource", "Cannot open file: %s",fname);
1862 if (!lenfile) delete [] fname;
1863 return;
1864 }
1865
1866 //set precision
1867 Int_t precision = gEnv->GetValue("Canvas.SavePrecision",7);
1868 out.precision(precision);
1869
1870 // Write macro header and date/time stamp
1871 TDatime t;
1873 Int_t topx,topy;
1874 UInt_t w, h;
1875 if (!fCanvasImp) {
1876 Error("SaveSource", "Cannot open TCanvas");
1877 return;
1878 }
1879 UInt_t editorWidth = fCanvasImp->GetWindowGeometry(topx,topy,w,h);
1880 w = UInt_t((fWindowWidth - editorWidth)/cx);
1881 h = UInt_t((fWindowHeight)/cx);
1882 topx = GetWindowTopX();
1883 topy = GetWindowTopY();
1884
1885 if (w == 0) {
1886 w = GetWw()+4; h = GetWh()+4;
1887 topx = 1; topy = 1;
1888 }
1889
1890 TString mname(fname);
1891 out << R"CODE(#ifdef __CLING__
1892#pragma cling optimize(0)
1893#endif
1894)CODE";
1895 Int_t p = mname.Last('.');
1896 Int_t s = mname.Last('/')+1;
1897
1898 // A named macro is generated only if the function name is valid. If not, the
1899 // macro is unnamed.
1900 TString first(mname(s,s+1));
1901 if (!first.IsDigit()) out <<"void " << mname(s,p-s) << "()" << std::endl;
1902
1903 out <<"{"<<std::endl;
1904 out <<"//=========Macro generated from canvas: "<<GetName()<<"/"<<GetTitle()<<std::endl;
1905 out <<"//========= ("<<t.AsString()<<") by ROOT version "<<gROOT->GetVersion()<<std::endl;
1906
1908 out <<std::endl<<" gStyle->SetCanvasPreferGL(kTRUE);"<<std::endl<<std::endl;
1909
1910 // Write canvas parameters (TDialogCanvas case)
1911 if (InheritsFrom(TDialogCanvas::Class())) {
1912 out<<" "<<ClassName()<<" *"<<cname<<" = new "<<ClassName()<<"("<<quote<<GetName()
1913 <<quote<<", "<<quote<<GetTitle()<<quote<<","<<w<<","<<h<<");"<<std::endl;
1914 } else {
1915 // Write canvas parameters (TCanvas case)
1916 out<<" TCanvas *"<<cname<<" = new TCanvas("<<quote<<GetName()<<quote<<", "<<quote<<GetTitle()
1917 <<quote;
1918 if (!HasMenuBar())
1919 out<<",-"<<topx<<","<<topy<<","<<w<<","<<h<<");"<<std::endl;
1920 else
1921 out<<","<<topx<<","<<topy<<","<<w<<","<<h<<");"<<std::endl;
1922 }
1923 // Write canvas options (in $TROOT or $TStyle)
1924 if (gStyle->GetOptFit()) {
1925 out<<" gStyle->SetOptFit(1);"<<std::endl;
1926 }
1927 if (!gStyle->GetOptStat()) {
1928 out<<" gStyle->SetOptStat(0);"<<std::endl;
1929 }
1930 if (!gStyle->GetOptTitle()) {
1931 out<<" gStyle->SetOptTitle(0);"<<std::endl;
1932 }
1933 if (gROOT->GetEditHistograms()) {
1934 out<<" gROOT->SetEditHistograms();"<<std::endl;
1935 }
1936 if (GetShowEventStatus()) {
1937 out<<" "<<GetName()<<"->ToggleEventStatus();"<<std::endl;
1938 }
1939 if (GetShowToolTips()) {
1940 out<<" "<<GetName()<<"->ToggleToolTips();"<<std::endl;
1941 }
1942 if (GetHighLightColor() != 5) {
1943 if (GetHighLightColor() > 228) {
1945 out<<" "<<GetName()<<"->SetHighLightColor(ci);" << std::endl;
1946 } else
1947 out<<" "<<GetName()<<"->SetHighLightColor("<<GetHighLightColor()<<");"<<std::endl;
1948 }
1949
1950 // Now recursively scan all pads of this canvas
1951 cd();
1952 if (invalid) SetName("c1");
1953 TPad::SavePrimitive(out,option);
1954 // Write canvas options related to pad editor
1955 out<<" "<<GetName()<<"->SetSelected("<<GetName()<<");"<<std::endl;
1956 if (GetShowToolBar()) {
1957 out<<" "<<GetName()<<"->ToggleToolBar();"<<std::endl;
1958 }
1959 if (invalid) SetName(" ");
1960
1961 out <<"}"<<std::endl;
1962 out.close();
1963 Info("SaveSource","C++ Macro file: %s has been generated", fname);
1964
1965 // reset bit TClass::kClassSaved for all classes
1966 next.Reset();
1967 while((cl = (TClass*)next())) {
1969 }
1970 if (!lenfile) delete [] fname;
1971}
1972
1973////////////////////////////////////////////////////////////////////////////////
1974/// Toggle batch mode. However, if the canvas is created without a window
1975/// then batch mode always stays set.
1976
1977void TCanvas::SetBatch(Bool_t batch)
1978{
1979 if (gROOT->IsBatch() || IsWeb())
1980 fBatch = kTRUE;
1981 else
1982 fBatch = batch;
1983}
1984
1985////////////////////////////////////////////////////////////////////////////////
1986/// Set Width and Height of canvas to ww and wh respectively. If ww and/or wh
1987/// are greater than the current canvas window a scroll bar is automatically
1988/// generated. Use this function to zoom in a canvas and navigate via
1989/// the scroll bars. The Width and Height in this method are different from those
1990/// given in the TCanvas constructors where these two dimension include the size
1991/// of the window decoration whereas they do not in this method.
1992
1994{
1995 if (fCanvasImp) {
1996 fCanvasImp->SetCanvasSize(ww, wh);
1997 fCw = ww;
1998 fCh = wh;
1999 TPad *padsav = (TPad*)gPad;
2000 cd();
2001 ResizePad();
2002 if (padsav) padsav->cd();
2003 }
2004}
2005
2006////////////////////////////////////////////////////////////////////////////////
2007/// Set cursor.
2008
2009void TCanvas::SetCursor(ECursor cursor)
2010{
2011 if (IsBatch()) return;
2012 gVirtualX->SetCursor(fCanvasID, cursor);
2013}
2014
2015////////////////////////////////////////////////////////////////////////////////
2016/// Set Double Buffer On/Off.
2017
2019{
2020 if (IsBatch()) return;
2021 fDoubleBuffer = mode;
2022 gVirtualX->SetDoubleBuffer(fCanvasID, mode);
2023
2024 // depending of the buffer mode set the drawing window to either
2025 // the canvas pixmap or to the canvas on-screen window
2026 if (fDoubleBuffer) {
2028 } else
2030}
2031
2032////////////////////////////////////////////////////////////////////////////////
2033/// Fix canvas aspect ratio to current value if fixed is true.
2034
2036{
2037 if (fixed) {
2038 if (!fFixedAspectRatio) {
2039 if (fCh != 0)
2041 else {
2042 Error("SetAspectRatio", "cannot fix aspect ratio, height of canvas is 0");
2043 return;
2044 }
2046 }
2047 } else {
2049 fAspectRatio = 0;
2050 }
2051}
2052
2053////////////////////////////////////////////////////////////////////////////////
2054/// If isfolder=kTRUE, the canvas can be browsed like a folder
2055/// by default a canvas is not browsable.
2056
2057void TCanvas::SetFolder(Bool_t isfolder)
2058{
2059 fgIsFolder = isfolder;
2060}
2061
2062////////////////////////////////////////////////////////////////////////////////
2063/// Set canvas name. In case `name` is an empty string, a default name is set.
2064
2065void TCanvas::SetName(const char *name)
2066{
2067 if (!name || !name[0]) {
2068 const char *defcanvas = gROOT->GetDefCanvasName();
2069 char *cdef;
2070 auto lc = (TList*)gROOT->GetListOfCanvases();
2071 if (lc->FindObject(defcanvas)) {
2072 cdef = Form("%s_n%d",defcanvas,lc->GetSize()+1);
2073 } else {
2074 cdef = Form("%s",defcanvas);
2075 }
2076 fName = cdef;
2077 } else {
2078 fName = name;
2079 }
2080 if (gPad && TestBit(kMustCleanup)) gPad->Modified();
2081}
2082
2083
2084////////////////////////////////////////////////////////////////////////////////
2085/// Function to resize a canvas so that the plot inside is shown in real aspect
2086/// ratio
2087///
2088/// \param[in] axis 1 for resizing horizontally (x-axis) in order to get real
2089/// aspect ratio, 2 for the resizing vertically (y-axis)
2090/// \return false if error is encountered, true otherwise
2091///
2092/// ~~~ {.cpp}
2093/// hpxpy->Draw();
2094/// c1->SetRealAspectRatio();
2095/// ~~~
2096///
2097/// - For defining the concept of real aspect ratio, it is assumed that x and y
2098/// axes are in same units, e.g. both in MeV or both in ns.
2099/// - You can resize either the width of the canvas or the height, but not both
2100/// at the same time
2101/// - Call this function AFTER drawing AND zooming (SetUserRange) your TGraph or
2102/// Histogram, otherwise it cannot infer your actual axes lengths
2103/// - This function ensures that the TFrame has a real aspect ratio, this does not
2104/// mean that the full pad (i.e. the canvas or png output) including margins has
2105/// exactly the same ratio
2106/// - This function does not work if the canvas is divided in several subpads
2107
2108bool TCanvas::SetRealAspectRatio(const Int_t axis)
2109{
2110 Update();
2111
2112 //Get how many pixels are occupied by the canvas
2113 Int_t npx = GetWw();
2114 Int_t npy = GetWh();
2115
2116 //Get x-y coordinates at the edges of the canvas (extrapolating outside the axes, NOT at the edges of the histogram)
2117 Double_t x1 = GetX1();
2118 Double_t y1 = GetY1();
2119 Double_t x2 = GetX2();
2120 Double_t y2 = GetY2();
2121
2122 //Get the length of extrapolated x and y axes
2123 Double_t xlength2 = x2 - x1;
2124 Double_t ylength2 = y2 - y1;
2125 Double_t ratio2 = xlength2/ylength2;
2126
2127 //Now get the number of pixels including the canvas borders
2128 Int_t bnpx = GetWindowWidth();
2129 Int_t bnpy = GetWindowHeight();
2130
2131 if (axis==1) {
2132 SetCanvasSize(TMath::Nint(npy*ratio2), npy);
2133 SetWindowSize((bnpx-npx)+TMath::Nint(npy*ratio2), bnpy);
2134 } else if (axis==2) {
2135 SetCanvasSize(npx, TMath::Nint(npx/ratio2));
2136 SetWindowSize(bnpx, (bnpy-npy)+TMath::Nint(npx/ratio2));
2137 } else {
2138 Error("SetRealAspectRatio", "axis value %d is neither 1 (resize along x-axis) nor 2 (resize along y-axis).",axis);
2139 return false;
2140 }
2141
2142 //Check now that resizing has worked
2143
2144 Update();
2145
2146 //Get how many pixels are occupied by the canvas
2147 npx = GetWw();
2148 npy = GetWh();
2149
2150 //Get x-y coordinates at the edges of the canvas (extrapolating outside the axes,
2151 //NOT at the edges of the histogram)
2152 x1 = GetX1();
2153 y1 = GetY1();
2154 x2 = GetX2();
2155 y2 = GetY2();
2156
2157 //Get the length of extrapolated x and y axes
2158 xlength2 = x2 - x1;
2159 ylength2 = y2 - y1;
2160 ratio2 = xlength2/ylength2;
2161
2162 //Check accuracy +/-1 pixel due to rounding
2163 if (abs(TMath::Nint(npy*ratio2) - npx)<2) {
2164 return true;
2165 } else {
2166 Error("SetRealAspectRatio", "Resizing failed.");
2167 return false;
2168 }
2169}
2170
2171
2172////////////////////////////////////////////////////////////////////////////////
2173/// Set selected canvas.
2174
2176{
2177 fSelected = obj;
2178 if (obj) obj->SetBit(kMustCleanup);
2179}
2180
2181////////////////////////////////////////////////////////////////////////////////
2182/// Set canvas title.
2183
2184void TCanvas::SetTitle(const char *title)
2185{
2186 fTitle = title;
2188}
2189
2190////////////////////////////////////////////////////////////////////////////////
2191/// Set canvas window position
2192
2194{
2195 if (fCanvasImp)
2197}
2198
2199////////////////////////////////////////////////////////////////////////////////
2200/// Set canvas window size
2201
2203{
2204 if (fBatch)
2205 SetCanvasSize((ww + fCw) / 2, (wh + fCh) / 2);
2206 else if (fCanvasImp)
2207 fCanvasImp->SetWindowSize(ww, wh);
2208}
2209
2210////////////////////////////////////////////////////////////////////////////////
2211/// Set the canvas scale in centimeters.
2212///
2213/// This information is used by PostScript to set the page size.
2214///
2215/// \param[in] xsize size of the canvas in centimeters along X
2216/// \param[in] ysize size of the canvas in centimeters along Y
2217///
2218/// if xsize and ysize are not equal to 0, then the scale factors will
2219/// be computed to keep the ratio ysize/xsize independently of the canvas
2220/// size (parts of the physical canvas will be unused).
2221///
2222/// if xsize = 0 and ysize is not zero, then xsize will be computed
2223/// to fit to the current canvas scale. If the canvas is resized,
2224/// a new value for xsize will be recomputed. In this case the aspect
2225/// ratio is not preserved.
2226///
2227/// if both xsize = 0 and ysize = 0, then the scaling is automatic.
2228/// the largest dimension will be allocated a size of 20 centimeters.
2229
2230void TCanvas::Size(Float_t xsize, Float_t ysize)
2231{
2232 fXsizeUser = xsize;
2233 fYsizeUser = ysize;
2234
2235 Resize();
2236}
2237
2238////////////////////////////////////////////////////////////////////////////////
2239/// Show canvas
2240
2241void TCanvas::Show()
2242{
2243 if (fCanvasImp)
2244 fCanvasImp->Show();
2245}
2246
2247////////////////////////////////////////////////////////////////////////////////
2248/// Stream a class object.
2249
2250void TCanvas::Streamer(TBuffer &b)
2251{
2252 UInt_t R__s, R__c;
2253 if (b.IsReading()) {
2254 Version_t v = b.ReadVersion(&R__s, &R__c);
2255 gPad = this;
2256 fCanvas = this;
2257 if (v>7) b.ClassBegin(TCanvas::IsA());
2258 if (v>7) b.ClassMember("TPad");
2259 TPad::Streamer(b);
2260 gPad = this;
2261 //restore the colors
2262 TObjArray *colors = (TObjArray*)fPrimitives->FindObject("ListOfColors");
2263 if (colors) {
2264 TIter next(colors);
2265 TColor *colold;
2266 while ((colold = (TColor*)next())) {
2267 if (colold) {
2268 Int_t cn = 0;
2269 if (colold) cn = colold->GetNumber();
2270 TColor *colcur = gROOT->GetColor(cn);
2271 if (colcur) {
2272 colcur->SetRGB(colold->GetRed(),colold->GetGreen(),colold->GetBlue());
2273 } else {
2274 colcur = new TColor(cn,colold->GetRed(),
2275 colold->GetGreen(),
2276 colold->GetBlue(),
2277 colold->GetName());
2278 if (!colcur) return;
2279 }
2280 }
2281 }
2282 //restore the palette if needed
2283 TObjArray *currentpalette = (TObjArray*)fPrimitives->FindObject("CurrentColorPalette");
2284 if (currentpalette) {
2285 TIter nextpal(currentpalette);
2286 Int_t n = currentpalette->GetEntries();
2287 TArrayI palcolors(n);
2288 TColor *col = 0;
2289 Int_t i = 0;
2290 while ((col = (TColor*)nextpal())) palcolors[i++] = col->GetNumber();
2291 gStyle->SetPalette(n,palcolors.GetArray());
2292 fPrimitives->Remove(currentpalette);
2293 delete currentpalette;
2294 }
2296 colors->Delete();
2297 delete colors;
2298 }
2299
2300 if (v>7) b.ClassMember("fDISPLAY","TString");
2301 fDISPLAY.Streamer(b);
2302 if (v>7) b.ClassMember("fDoubleBuffer", "Int_t");
2303 b >> fDoubleBuffer;
2304 if (v>7) b.ClassMember("fRetained", "Bool_t");
2305 b >> fRetained;
2306 if (v>7) b.ClassMember("fXsizeUser", "Size_t");
2307 b >> fXsizeUser;
2308 if (v>7) b.ClassMember("fYsizeUser", "Size_t");
2309 b >> fYsizeUser;
2310 if (v>7) b.ClassMember("fXsizeReal", "Size_t");
2311 b >> fXsizeReal;
2312 if (v>7) b.ClassMember("fYsizeReal", "Size_t");
2313 b >> fYsizeReal;
2314 fCanvasID = -1;
2315 if (v>7) b.ClassMember("fWindowTopX", "Int_t");
2316 b >> fWindowTopX;
2317 if (v>7) b.ClassMember("fWindowTopY", "Int_t");
2318 b >> fWindowTopY;
2319 if (v > 2) {
2320 if (v>7) b.ClassMember("fWindowWidth", "UInt_t");
2321 b >> fWindowWidth;
2322 if (v>7) b.ClassMember("fWindowHeight", "UInt_t");
2323 b >> fWindowHeight;
2324 }
2325 if (v>7) b.ClassMember("fCw", "UInt_t");
2326 b >> fCw;
2327 if (v>7) b.ClassMember("fCh", "UInt_t");
2328 b >> fCh;
2329 if (v <= 2) {
2330 fWindowWidth = fCw;
2332 }
2333 if (v>7) b.ClassMember("fCatt", "TAttCanvas");
2334 fCatt.Streamer(b);
2335 Bool_t dummy;
2336 if (v>7) b.ClassMember("kMoveOpaque", "Bool_t");
2337 b >> dummy; if (dummy) MoveOpaque(1);
2338 if (v>7) b.ClassMember("kResizeOpaque", "Bool_t");
2339 b >> dummy; if (dummy) ResizeOpaque(1);
2340 if (v>7) b.ClassMember("fHighLightColor", "Color_t");
2341 b >> fHighLightColor;
2342 if (v>7) b.ClassMember("fBatch", "Bool_t");
2343 b >> dummy; //was fBatch
2344 if (v < 2) return;
2345 if (v>7) b.ClassMember("kShowEventStatus", "Bool_t");
2346 b >> dummy; if (dummy) SetBit(kShowEventStatus);
2347
2348 if (v > 3) {
2349 if (v>7) b.ClassMember("kAutoExec", "Bool_t");
2350 b >> dummy; if (dummy) SetBit(kAutoExec);
2351 }
2352 if (v>7) b.ClassMember("kMenuBar", "Bool_t");
2353 b >> dummy; if (dummy) SetBit(kMenuBar);
2354 fBatch = gROOT->IsBatch();
2355 if (v>7) b.ClassEnd(TCanvas::IsA());
2356 b.CheckByteCount(R__s, R__c, TCanvas::IsA());
2357 } else {
2358 //save list of colors
2359 //we must protect the case when two or more canvases are saved
2360 //in the same buffer. If the list of colors has already been saved
2361 //in the buffer, do not add the list of colors to the list of primitives.
2362 TObjArray *colors = 0;
2363 TObjArray *CurrentColorPalette = 0;
2364 if (TColor::DefinedColors()) {
2365 if (!b.CheckObject(gROOT->GetListOfColors(),TObjArray::Class())) {
2366 colors = (TObjArray*)gROOT->GetListOfColors();
2368 }
2369 //save the current palette
2371 Int_t palsize = pal.GetSize();
2372 CurrentColorPalette = new TObjArray();
2373 CurrentColorPalette->SetName("CurrentColorPalette");
2374 for (Int_t i=0; i<palsize; i++) CurrentColorPalette->Add(gROOT->GetColor(pal[i]));
2375 fPrimitives->Add(CurrentColorPalette);
2376 }
2377
2378 R__c = b.WriteVersion(TCanvas::IsA(), kTRUE);
2379 b.ClassBegin(TCanvas::IsA());
2380 b.ClassMember("TPad");
2381 TPad::Streamer(b);
2383 if (CurrentColorPalette) { fPrimitives->Remove(CurrentColorPalette); delete CurrentColorPalette; }
2384 b.ClassMember("fDISPLAY","TString");
2385 fDISPLAY.Streamer(b);
2386 b.ClassMember("fDoubleBuffer", "Int_t");
2387 b << fDoubleBuffer;
2388 b.ClassMember("fRetained", "Bool_t");
2389 b << fRetained;
2390 b.ClassMember("fXsizeUser", "Size_t");
2391 b << fXsizeUser;
2392 b.ClassMember("fYsizeUser", "Size_t");
2393 b << fYsizeUser;
2394 b.ClassMember("fXsizeReal", "Size_t");
2395 b << fXsizeReal;
2396 b.ClassMember("fYsizeReal", "Size_t");
2397 b << fYsizeReal;
2399 Int_t topx = fWindowTopX, topy = fWindowTopY;
2400 UInt_t editorWidth = 0;
2401 if(fCanvasImp) editorWidth = fCanvasImp->GetWindowGeometry(topx,topy,w,h);
2402 b.ClassMember("fWindowTopX", "Int_t");
2403 b << topx;
2404 b.ClassMember("fWindowTopY", "Int_t");
2405 b << topy;
2406 b.ClassMember("fWindowWidth", "UInt_t");
2407 b << (UInt_t)(w-editorWidth);
2408 b.ClassMember("fWindowHeight", "UInt_t");
2409 b << h;
2410 b.ClassMember("fCw", "UInt_t");
2411 b << fCw;
2412 b.ClassMember("fCh", "UInt_t");
2413 b << fCh;
2414 b.ClassMember("fCatt", "TAttCanvas");
2415 fCatt.Streamer(b);
2416 b.ClassMember("kMoveOpaque", "Bool_t");
2417 b << TestBit(kMoveOpaque); //please remove in ROOT version 6
2418 b.ClassMember("kResizeOpaque", "Bool_t");
2419 b << TestBit(kResizeOpaque); //please remove in ROOT version 6
2420 b.ClassMember("fHighLightColor", "Color_t");
2421 b << fHighLightColor;
2422 b.ClassMember("fBatch", "Bool_t");
2423 b << fBatch; //please remove in ROOT version 6
2424 b.ClassMember("kShowEventStatus", "Bool_t");
2425 b << TestBit(kShowEventStatus); //please remove in ROOT version 6
2426 b.ClassMember("kAutoExec", "Bool_t");
2427 b << TestBit(kAutoExec); //please remove in ROOT version 6
2428 b.ClassMember("kMenuBar", "Bool_t");
2429 b << TestBit(kMenuBar); //please remove in ROOT version 6
2430 b.ClassEnd(TCanvas::IsA());
2431 b.SetByteCount(R__c, kTRUE);
2432 }
2433}
2434
2435////////////////////////////////////////////////////////////////////////////////
2436/// Toggle pad auto execution of list of TExecs.
2437
2439{
2440 Bool_t autoExec = TestBit(kAutoExec);
2441 SetBit(kAutoExec,!autoExec);
2442}
2443
2444////////////////////////////////////////////////////////////////////////////////
2445/// Toggle event statusbar.
2446
2448{
2449 Bool_t showEventStatus = !TestBit(kShowEventStatus);
2450 SetBit(kShowEventStatus,showEventStatus);
2451
2452 if (fCanvasImp) fCanvasImp->ShowStatusBar(showEventStatus);
2453}
2454
2455////////////////////////////////////////////////////////////////////////////////
2456/// Toggle toolbar.
2457
2459{
2460 Bool_t showToolBar = !TestBit(kShowToolBar);
2461 SetBit(kShowToolBar,showToolBar);
2462
2463 if (fCanvasImp) fCanvasImp->ShowToolBar(showToolBar);
2464}
2465
2466////////////////////////////////////////////////////////////////////////////////
2467/// Toggle editor.
2468
2470{
2471 Bool_t showEditor = !TestBit(kShowEditor);
2472 SetBit(kShowEditor,showEditor);
2473
2474 if (fCanvasImp) fCanvasImp->ShowEditor(showEditor);
2475}
2476
2477////////////////////////////////////////////////////////////////////////////////
2478/// Toggle tooltip display.
2479
2481{
2482 Bool_t showToolTips = !TestBit(kShowToolTips);
2483 SetBit(kShowToolTips, showToolTips);
2484
2485 if (fCanvasImp) fCanvasImp->ShowToolTips(showToolTips);
2486}
2487
2488
2489////////////////////////////////////////////////////////////////////////////////
2490/// Static function returning "true" if transparency is supported.
2491
2493{
2494 return gPad && (gVirtualX->InheritsFrom("TGQuartz") ||
2495 gPad->GetGLDevice() != -1);
2496}
2497
2498extern "C" void ROOT_TCanvas_Update(void* TheCanvas) {
2499 static_cast<TCanvas*>(TheCanvas)->Update();
2500}
2501
2502////////////////////////////////////////////////////////////////////////////////
2503/// Update canvas pad buffers.
2504
2505void TCanvas::Update()
2506{
2507 if (fUpdating) return;
2508
2509 if (fPixmapID == -1) return;
2510
2511 static const union CastFromFuncToVoidPtr_t {
2512 CastFromFuncToVoidPtr_t(): fFuncPtr(ROOT_TCanvas_Update) {}
2513 void (*fFuncPtr)(void*);
2514 void* fVoidPtr;
2515 } castFromFuncToVoidPtr;
2516
2517 if (gThreadXAR) {
2518 void *arr[3];
2519 arr[1] = this;
2520 arr[2] = castFromFuncToVoidPtr.fVoidPtr;
2521 if ((*gThreadXAR)("CUPD", 3, arr, 0)) return;
2522 }
2523
2524 if (!fCanvasImp) return;
2525
2526 if (!gVirtualX->IsCmdThread()) {
2527 // Why do we have this (which uses the interpreter to funnel the Update()
2528 // through the main thread) when the gThreadXAR mechanism does seemingly
2529 // the same?
2530 gInterpreter->Execute(this, IsA(), "Update", "");
2531 return;
2532 }
2533
2535
2536 fUpdating = kTRUE;
2537
2538 if (!fCanvasImp->PerformUpdate()) {
2539
2540 if (!IsBatch()) FeedbackMode(kFALSE); // Goto double buffer mode
2541
2542 if (!UseGL() || fGLDevice == -1) PaintModified(); // Repaint all modified pad's
2543
2544 Flush(); // Copy all pad pixmaps to the screen
2545
2547 }
2548
2549 fUpdating = kFALSE;
2550}
2551
2552////////////////////////////////////////////////////////////////////////////////
2553/// Used by friend class TCanvasImp.
2554
2556{
2557 fCanvasID = 0;
2558 fContextMenu = 0;
2559}
2560
2561////////////////////////////////////////////////////////////////////////////////
2562/// Check whether this canvas is to be drawn in grayscale mode.
2563
2565{
2566 return TestBit(kIsGrayscale);
2567}
2568
2569////////////////////////////////////////////////////////////////////////////////
2570/// Set whether this canvas should be painted in grayscale, and re-paint
2571/// it if necessary.
2572
2573void TCanvas::SetGrayscale(Bool_t set /*= kTRUE*/)
2574{
2575 if (IsGrayscale() == set) return;
2576 SetBit(kIsGrayscale, set);
2577 Paint(); // update canvas and all sub-pads, unconditionally!
2578}
2579
2580////////////////////////////////////////////////////////////////////////////////
2581/// Probably, TPadPainter must be placed in a separate ROOT module -
2582/// "padpainter" (the same as "histpainter"). But now, it's directly in a
2583/// gpad dir, so, in case of default painter, no *.so should be loaded,
2584/// no need in plugin managers.
2585/// May change in future.
2586
2588{
2589 //Even for batch mode painter is still required, just to delegate
2590 //some calls to batch "virtual X".
2591 if (!UseGL() || fBatch) {
2592 fPainter = 0;
2594 if (!fPainter) fPainter = new TPadPainter; // Do not need plugin manager for this!
2595 } else {
2597 if (!fPainter) {
2598 Error("CreatePainter", "GL Painter creation failed! Will use default!");
2599 fPainter = new TPadPainter;
2600 fUseGL = kFALSE;
2601 }
2602 }
2603}
2604
2605////////////////////////////////////////////////////////////////////////////////
2606/// Access and (probably) creation of pad painter.
2607
2609{
2610 if (!fPainter) CreatePainter();
2611 return fPainter;
2612}
2613
2614
2615////////////////////////////////////////////////////////////////////////////////
2616///assert on IsBatch() == false?
2617
2619{
2620 if (fGLDevice != -1) {
2621 //fPainter has a font manager.
2622 //Font manager will delete textures.
2623 //If context is wrong (we can have several canvases) -
2624 //wrong texture will be deleted, damaging some of our fonts.
2625 gGLManager->MakeCurrent(fGLDevice);
2626 }
2627
2629
2630 if (fGLDevice != -1) {
2631 gGLManager->DeleteGLContext(fGLDevice);//?
2632 fGLDevice = -1;
2633 }
2634}
EEventType
Definition Buttons.h:15
@ kButton1ShiftMotion
Definition Buttons.h:18
@ kMouseMotion
Definition Buttons.h:23
@ kWheelUp
Definition Buttons.h:18
@ kButton3Up
Definition Buttons.h:19
@ kButton2Motion
Definition Buttons.h:20
@ kButton3Motion
Definition Buttons.h:20
@ kButton3Down
Definition Buttons.h:17
@ kButton2Down
Definition Buttons.h:17
@ kKeyPress
Definition Buttons.h:20
@ kButton2Double
Definition Buttons.h:24
@ kArrowKeyRelease
Definition Buttons.h:21
@ kButton1Double
Definition Buttons.h:24
@ kButton3Double
Definition Buttons.h:24
@ kButton1Shift
Definition Buttons.h:18
@ kButton1Motion
Definition Buttons.h:20
@ kButton1Up
Definition Buttons.h:19
@ kWheelDown
Definition Buttons.h:18
@ kArrowKeyPress
Definition Buttons.h:21
@ kButton2Up
Definition Buttons.h:19
@ kMouseLeave
Definition Buttons.h:23
@ kButton1Down
Definition Buttons.h:17
@ kMouseEnter
Definition Buttons.h:23
typedef void(GLAPIENTRYP _GLUfuncptr)(void)
ECursor
Definition GuiTypes.h:372
@ kCross
Definition GuiTypes.h:374
#define SafeDelete(p)
Definition RConfig.hxx:537
#define b(i)
Definition RSha256.hxx:100
#define c(i)
Definition RSha256.hxx:101
#define h(i)
Definition RSha256.hxx:106
static const double x2[5]
static const double x1[5]
int Int_t
Definition RtypesCore.h:45
float Size_t
Definition RtypesCore.h:96
long Longptr_t
Definition RtypesCore.h:82
short Version_t
Definition RtypesCore.h:65
unsigned int UInt_t
Definition RtypesCore.h:46
const Bool_t kFALSE
Definition RtypesCore.h:101
bool Bool_t
Definition RtypesCore.h:63
double Double_t
Definition RtypesCore.h:59
float Float_t
Definition RtypesCore.h:57
const Bool_t kTRUE
Definition RtypesCore.h:100
const char Option_t
Definition RtypesCore.h:66
@ kRed
Definition Rtypes.h:66
R__EXTERN TApplication * gApplication
void ROOT_TCanvas_Update(void *TheCanvas)
Definition TCanvas.cxx:2495
class TCanvasInit gCanvasInit
const Size_t kDefaultCanvasSize
Definition TCanvas.cxx:61
R__EXTERN TEnv * gEnv
Definition TEnv.h:170
char name[80]
Definition TGX11.cxx:110
R__EXTERN TGuiFactory * gBatchGuiFactory
Definition TGuiFactory.h:67
R__EXTERN TGuiFactory * gGuiFactory
Definition TGuiFactory.h:66
#define gInterpreter
#define ClassImpQ(name)
Definition TQObject.h:282
Int_t gDebug
Definition TROOT.cxx:592
R__EXTERN TVirtualMutex * gROOTMutex
Definition TROOT.h:63
#define gROOT
Definition TROOT.h:404
char * Form(const char *fmt,...)
char * StrDup(const char *str)
Duplicate the string str.
Definition TString.cxx:2515
R__EXTERN TStyle * gStyle
Definition TStyle.h:413
#define gGLManager
Definition TVirtualGL.h:162
#define R__LOCKGUARD(mutex)
R__EXTERN TVirtualPS * gVirtualPS
Definition TVirtualPS.h:81
#define gPad
R__EXTERN Int_t(* gThreadXAR)(const char *xact, Int_t nb, void **ar, Int_t *iret)
#define gVirtualX
Definition TVirtualX.h:338
Color * colors
Definition X3DBuffer.c:21
#define snprintf
Definition civetweb.c:1540
void InitializeGraphics()
Initialize the graphics environment.
static void CreateApplication()
Static function used to create a default application environment.
static void NeedGraphicsLibs()
Static method.
Array of integers (32 bits per element).
Definition TArrayI.h:27
Int_t GetSize() const
Definition TArray.h:47
Fill Area Attributes class.
Definition TAttFill.h:19
virtual Color_t GetFillColor() const
Return the fill area color.
Definition TAttFill.h:30
void Copy(TAttFill &attfill) const
Copy this fill attributes to a new TAttFill.
Definition TAttFill.cxx:204
virtual void SetFillColor(Color_t fcolor)
Set the fill area color.
Definition TAttFill.h:37
Line Attributes class.
Definition TAttLine.h:18
void Copy(TAttLine &attline) const
Copy this line attributes to a new TAttLine.
Definition TAttLine.cxx:175
Manages default Pad attributes.
Definition TAttPad.h:19
virtual void SetBottomMargin(Float_t bottommargin)
Set Pad bottom margin in fraction of the pad height.
Definition TAttPad.cxx:99
virtual void SetLeftMargin(Float_t leftmargin)
Set Pad left margin in fraction of the pad width.
Definition TAttPad.cxx:109
virtual void Copy(TAttPad &attpad) const
copy function
Definition TAttPad.cxx:44
virtual void SetRightMargin(Float_t rightmargin)
Set Pad right margin in fraction of the pad width.
Definition TAttPad.cxx:119
virtual void SetTopMargin(Float_t topmargin)
Set Pad top margin in fraction of the pad height.
Definition TAttPad.cxx:129
Class to manage histogram axis.
Definition TAxis.h:30
virtual Bool_t GetTimeDisplay() const
Definition TAxis.h:126
virtual const char * GetTimeFormat() const
Definition TAxis.h:127
Create a Box.
Definition TBox.h:22
Using a TBrowser one can browse all ROOT objects.
Definition TBrowser.h:37
Buffer base class used for serializing objects.
Definition TBuffer.h:43
virtual void Show()
Definition TCanvasImp.h:66
virtual void Iconify()
Definition TCanvasImp.h:59
virtual Int_t InitWindow()
Definition TCanvasImp.h:60
virtual void Close()
Definition TCanvasImp.h:56
virtual void SetWindowPosition(Int_t x, Int_t y)
Definition TCanvasImp.h:90
virtual void SetWindowTitle(const char *newTitle)
Definition TCanvasImp.h:92
virtual Bool_t IsWeb() const
Definition TCanvasImp.h:45
virtual TVirtualPadPainter * CreatePadPainter()
Definition TCanvasImp.h:47
virtual void ShowToolBar(Bool_t show=kTRUE)
Definition TCanvasImp.h:100
virtual void RaiseWindow()
Definition TCanvasImp.h:96
virtual void ShowMenuBar(Bool_t show=kTRUE)
Definition TCanvasImp.h:94
virtual void ShowEditor(Bool_t show=kTRUE)
Definition TCanvasImp.h:99
virtual Bool_t PerformUpdate()
Definition TCanvasImp.h:46
virtual void SetWindowSize(UInt_t w, UInt_t h)
Definition TCanvasImp.h:91
virtual void ShowStatusBar(Bool_t show=kTRUE)
Definition TCanvasImp.h:95
virtual void ShowToolTips(Bool_t show=kTRUE)
Definition TCanvasImp.h:101
virtual void SetStatusText(const char *text=0, Int_t partidx=0)
Definition TCanvasImp.h:89
virtual UInt_t GetWindowGeometry(Int_t &x, Int_t &y, UInt_t &w, UInt_t &h)
Definition TCanvasImp.h:87
virtual void ForceUpdate()
Definition TCanvasImp.h:57
virtual void SetCanvasSize(UInt_t w, UInt_t h)
Definition TCanvasImp.h:93
The Canvas class.
Definition TCanvas.h:23
void Init()
Initialize the TCanvas members. Called by all constructors.
Definition TCanvas.cxx:525
UInt_t fCw
Width of the canvas along X (pixels)
Definition TCanvas.h:42
void EmbedInto(Int_t winid, Int_t ww, Int_t wh)
Embedded a canvas into a TRootEmbeddedCanvas.
Definition TCanvas.cxx:1037
void SetWindowSize(UInt_t ww, UInt_t wh)
Set canvas window size.
Definition TCanvas.cxx:2199
static void SetFolder(Bool_t isfolder=kTRUE)
If isfolder=kTRUE, the canvas can be browsed like a folder by default a canvas is not browsable.
Definition TCanvas.cxx:2054
void Browse(TBrowser *b) override
Browse.
Definition TCanvas.cxx:672
UInt_t GetWindowHeight() const
Definition TCanvas.h:160
virtual void EditorBar()
Get editor bar.
Definition TCanvas.cxx:1028
static TCanvas * MakeDefCanvas()
Static function to build a default canvas.
Definition TCanvas.cxx:1504
void EnterLeave(TPad *prevSelPad, TObject *prevSelObj)
Generate kMouseEnter and kMouseLeave events depending on the previously selected object and the curre...
Definition TCanvas.cxx:1063
Size_t fYsizeReal
Current size of canvas along Y in CM.
Definition TCanvas.h:35
void Constructor()
Canvas default constructor.
Definition TCanvas.cxx:181
virtual void ToggleAutoExec()
Toggle pad auto execution of list of TExecs.
Definition TCanvas.cxx:2435
TCanvas(const TCanvas &canvas)=delete
Int_t fWindowTopX
Top X position of window (in pixels)
Definition TCanvas.h:38
void Draw(Option_t *option="") override
Draw a canvas.
Definition TCanvas.cxx:843
void SetDoubleBuffer(Int_t mode=1) override
Set Double Buffer On/Off.
Definition TCanvas.cxx:2015
virtual void ToggleToolTips()
Toggle tooltip display.
Definition TCanvas.cxx:2477
void Clear(Option_t *option="") override
Remove all primitives from the canvas.
Definition TCanvas.cxx:725
void UseCurrentStyle() override
Force a copy of current style for all objects in canvas.
Definition TCanvas.cxx:1174
void Iconify()
Iconify canvas.
Definition TCanvas.cxx:1467
Int_t GetWindowTopX()
Returns current top x position of window on screen.
Definition TCanvas.cxx:1199
Bool_t IsWeb() const
Is web canvas.
Definition TCanvas.cxx:1484
virtual void ToggleEventStatus()
Toggle event statusbar.
Definition TCanvas.cxx:2444
void Destructor()
Actual canvas destructor.
Definition TCanvas.cxx:682
void DeleteCanvasPainter()
assert on IsBatch() == false?
Definition TCanvas.cxx:2615
TPad * fPadSave
! Pointer to saved pad in HandleInput
Definition TCanvas.h:55
static Bool_t SupportAlpha()
Static function returning "true" if transparency is supported.
Definition TCanvas.cxx:2489
Bool_t fBatch
! True when in batchmode
Definition TCanvas.h:58
Bool_t fUseGL
! True when rendering is with GL
Definition TCanvas.h:61
Int_t fEventX
! Last X mouse position in canvas
Definition TCanvas.h:45
Bool_t IsBatch() const override
Definition TCanvas.h:169
TObject * DrawClone(Option_t *option="") const override
Draw a clone of this canvas A new canvas is created that is a clone of this canvas.
Definition TCanvas.cxx:890
Size_t fXsizeReal
Current size of canvas along X in CM.
Definition TCanvas.h:34
Bool_t HasMenuBar() const
Definition TCanvas.h:166
TVirtualPadPainter * GetCanvasPainter()
Access and (probably) creation of pad painter.
Definition TCanvas.cxx:2605
virtual void HighlightConnect(const char *slot)
This is "simplification" for function TCanvas::Connect with Highlighted signal for specific slot.
Definition TCanvas.cxx:1631
TPad * Pick(Int_t px, Int_t py, TObjLink *&pickobj) override
Search for an object at pixel position px,py.
Definition TCanvas.h:180
void Close(Option_t *option="") override
Close canvas.
Definition TCanvas.cxx:776
void SetFixedAspectRatio(Bool_t fixed=kTRUE) override
Fix canvas aspect ratio to current value if fixed is true.
Definition TCanvas.cxx:2032
virtual void Resize(Option_t *option="")
Recompute canvas parameters following a X11 Resize.
Definition TCanvas.cxx:1668
Color_t GetHighLightColor() const override
Definition TCanvas.h:136
Bool_t GetShowToolBar() const
Definition TCanvas.h:147
void DrawEventStatus(Int_t event, Int_t x, Int_t y, TObject *selected)
Report name and title of primitive below the cursor.
Definition TCanvas.cxx:962
Bool_t IsFolder() const override
Is folder ?
Definition TCanvas.cxx:1476
UInt_t fWindowWidth
Width of window (including borders, etc.)
Definition TCanvas.h:40
TVirtualPadPainter * fPainter
! Canvas (pad) painter.
Definition TCanvas.h:64
void CopyPixmaps() override
Copy the canvas pixmap of the pad to the canvas.
Definition TCanvas.cxx:825
Bool_t IsGrayscale()
Check whether this canvas is to be drawn in grayscale mode.
Definition TCanvas.cxx:2561
TPad * fClickSelectedPad
! Pad containing currently click-selected object
Definition TCanvas.h:54
Bool_t fUpdating
! True when Updating the canvas
Definition TCanvas.h:59
void SaveSource(const char *filename="", Option_t *option="")
Save primitives in this canvas as a C++ macro file.
Definition TCanvas.cxx:1826
Color_t fHighLightColor
Highlight color of active pad.
Definition TCanvas.h:36
virtual void Size(Float_t xsizeuser=0, Float_t ysizeuser=0)
Set the canvas scale in centimeters.
Definition TCanvas.cxx:2227
virtual void ProcessedEvent(Int_t event, Int_t x, Int_t y, TObject *selected)
Emit ProcessedEvent() signal.
Definition TCanvas.cxx:1653
virtual void HandleInput(EEventType button, Int_t x, Int_t y)
Handle Input Events.
Definition TCanvas.cxx:1223
Size_t fXsizeUser
User specified size of canvas along X in CM.
Definition TCanvas.h:32
Int_t fEventY
! Last Y mouse position in canvas
Definition TCanvas.h:46
TVirtualPad * cd(Int_t subpadnumber=0) override
Set current canvas & pad.
Definition TCanvas.cxx:706
UInt_t fWindowHeight
Height of window (including menubar, borders, etc.)
Definition TCanvas.h:41
Int_t GetWindowTopY()
Returns current top y position of window on screen.
Definition TCanvas.cxx:1210
TObject * fClickSelected
! Currently click-selected object
Definition TCanvas.h:49
void SetCanvasSize(UInt_t ww, UInt_t wh) override
Set Width and Height of canvas to ww and wh respectively.
Definition TCanvas.cxx:1990
void Show()
Show canvas.
Definition TCanvas.cxx:2238
TPad * fSelectedPad
! Pad containing currently selected object
Definition TCanvas.h:53
virtual void Selected(TVirtualPad *pad, TObject *obj, Int_t event)
Emit Selected() signal.
Definition TCanvas.cxx:1639
Int_t fSelectedX
! X of selected object
Definition TCanvas.h:50
virtual void ToggleEditor()
Toggle editor.
Definition TCanvas.cxx:2466
TVirtualPad * GetSelectedPad() const override
Definition TCanvas.h:144
virtual void Picked(TPad *selpad, TObject *selected, Int_t event)
Emit Picked() signal.
Definition TCanvas.cxx:1593
TObject * fSelected
! Currently selected object
Definition TCanvas.h:48
void SavePrimitive(std::ostream &out, Option_t *option="") override
Save primitives in this canvas in C++ macro file with GUI.
Definition TCanvas.cxx:1781
void SetCursor(ECursor cursor) override
Set cursor.
Definition TCanvas.cxx:2006
Bool_t GetShowToolTips() const
Definition TCanvas.h:149
Int_t fCanvasID
! Canvas identifier
Definition TCanvas.h:47
void SetGrayscale(Bool_t set=kTRUE)
Set whether this canvas should be painted in grayscale, and re-paint it if necessary.
Definition TCanvas.cxx:2570
void SetTitle(const char *title="") override
Set canvas title.
Definition TCanvas.cxx:2181
UInt_t fCh
Height of the canvas along Y (pixels)
Definition TCanvas.h:43
TContextMenu * fContextMenu
! Context menu pointer
Definition TCanvas.h:57
TAttCanvas fCatt
Canvas attributes.
Definition TCanvas.h:30
void SetName(const char *name="") override
Set canvas name. In case name is an empty string, a default name is set.
Definition TCanvas.cxx:2062
UInt_t GetWindowWidth() const
Definition TCanvas.h:159
Bool_t fRetained
Retain structure flag.
Definition TCanvas.h:60
void DisconnectWidget()
Used by friend class TCanvasImp.
Definition TCanvas.cxx:2552
void FeedbackMode(Bool_t set)
Turn rubberband feedback mode on or off.
Definition TCanvas.cxx:1118
void ls(Option_t *option="") const override
List all pads.
Definition TCanvas.cxx:1492
void RaiseWindow()
Raise canvas window.
Definition TCanvas.cxx:1747
void Build()
Build a canvas. Called by all constructors.
Definition TCanvas.cxx:574
Int_t fWindowTopY
Top Y position of window (in pixels)
Definition TCanvas.h:39
void Paint(Option_t *option="") override
Paint canvas.
Definition TCanvas.cxx:1543
void Update() override
Update canvas pad buffers.
Definition TCanvas.cxx:2502
void ExecuteEvent(Int_t event, Int_t px, Int_t py) override
Execute action corresponding to one event.
Definition TCanvas.cxx:1100
void RunAutoExec()
Execute the list of TExecs in the current pad.
Definition TCanvas.cxx:1770
virtual void Cleared(TVirtualPad *pad)
Emit pad Cleared signal.
Definition TCanvas.cxx:758
UInt_t GetWw() const override
Definition TCanvas.h:161
TCanvasImp * fCanvasImp
! Window system specific canvas implementation
Definition TCanvas.h:56
UInt_t GetWh() const override
Definition TCanvas.h:162
virtual void Highlighted(TVirtualPad *pad, TObject *obj, Int_t x, Int_t y)
Emit Highlighted() signal.
Definition TCanvas.cxx:1612
void Flush()
Flush canvas buffers.
Definition TCanvas.cxx:1132
Size_t fYsizeUser
User specified size of canvas along Y in CM.
Definition TCanvas.h:33
Int_t fDoubleBuffer
Double buffer flag (0=off, 1=on)
Definition TCanvas.h:37
void ForceUpdate()
Force canvas update.
Definition TCanvas.cxx:1166
void CreatePainter()
Probably, TPadPainter must be placed in a separate ROOT module - "padpainter" (the same as "histpaint...
Definition TCanvas.cxx:2584
void SetSelected(TObject *obj) override
Set selected canvas.
Definition TCanvas.cxx:2172
void MoveOpaque(Int_t set=1)
Set option to move objects/pads in a canvas.
Definition TCanvas.cxx:1535
virtual ~TCanvas()
Canvas destructor.
Definition TCanvas.cxx:664
static Bool_t fgIsFolder
Indicates if canvas can be browsed as a folder.
Definition TCanvas.h:66
void Closed() override
Emit Closed signal.
Definition TCanvas.cxx:766
void SetWindowPosition(Int_t x, Int_t y)
Set canvas window position.
Definition TCanvas.cxx:2190
TString fDISPLAY
Name of destination screen.
Definition TCanvas.h:31
bool SetRealAspectRatio(const Int_t axis=1)
Function to resize a canvas so that the plot inside is shown in real aspect ratio.
Definition TCanvas.cxx:2105
Int_t fEvent
! Type of current or last handled event
Definition TCanvas.h:44
Bool_t GetShowEventStatus() const
Definition TCanvas.h:146
TString fSelectedOpt
! Drawing option of selected object
Definition TCanvas.h:52
Int_t fSelectedY
! Y of selected object
Definition TCanvas.h:51
Bool_t fDrawn
! Set to True when the Draw method is called
Definition TCanvas.h:62
void SetBatch(Bool_t batch=kTRUE) override
Toggle batch mode.
Definition TCanvas.cxx:1974
Bool_t UseGL() const
Definition TCanvas.h:223
@ kResizeOpaque
Definition TCanvas.h:93
@ kShowToolTips
Definition TCanvas.h:95
@ kShowToolBar
Definition TCanvas.h:90
@ kMoveOpaque
Definition TCanvas.h:92
@ kIsGrayscale
Definition TCanvas.h:94
@ kShowEventStatus
Definition TCanvas.h:87
@ kAutoExec
Definition TCanvas.h:88
@ kMenuBar
Definition TCanvas.h:89
@ kShowEditor
Definition TCanvas.h:91
void ResizeOpaque(Int_t set=1)
Set option to resize objects/pads in a canvas.
Definition TCanvas.cxx:1762
virtual void ToggleToolBar()
Toggle toolbar.
Definition TCanvas.cxx:2455
virtual TObject * DrawClonePad()
Draw a clone of this canvas into the current pad In an interactive session, select the destination/cu...
Definition TCanvas.cxx:907
TClass instances represent classes, structs and namespaces in the ROOT type system.
Definition TClass.h:80
@ kRealNew
Definition TClass.h:107
static ENewType IsCallingNew()
Static method returning the defConstructor flag passed to TClass::New().
Definition TClass.cxx:5888
@ kClassSaved
Definition TClass.h:94
void SetName(const char *name)
void Browse(TBrowser *b)
Browse this collection (called by TBrowser).
The color creation and management class.
Definition TColor.h:19
virtual void SetRGB(Float_t r, Float_t g, Float_t b)
Initialize this color and its associated colors.
Definition TColor.cxx:1758
static const TArrayI & GetPalette()
Static function returning the current active palette.
Definition TColor.cxx:1461
Float_t GetRed() const
Definition TColor.h:58
static Int_t GetColor(const char *hexcolor)
Static method returning color number for color specified by hex color string of form: "#rrggbb",...
Definition TColor.cxx:1822
static void SaveColor(std::ostream &out, Int_t ci)
Save a color with index > 228 as a C++ statement(s) on output stream out.
Definition TColor.cxx:2174
Int_t GetNumber() const
Definition TColor.h:56
Float_t GetBlue() const
Definition TColor.h:60
Float_t GetGreen() const
Definition TColor.h:59
static Bool_t DefinedColors()
Static function returning kTRUE if some new colors have been defined after initialisation or since th...
Definition TColor.cxx:1479
This class provides an interface to context sensitive popup menus.
virtual void Popup(Int_t x, Int_t y, TObject *obj, TVirtualPad *c=nullptr, TVirtualPad *p=nullptr)
Popup context menu at given location in canvas c and pad p for selected object.
This class stores the date and time with a precision of one second in an unsigned 32 bit word (950130...
Definition TDatime.h:37
const char * AsSQLString() const
Return the date & time in SQL compatible string format, like: 1997-01-15 20:16:28.
Definition TDatime.cxx:152
UInt_t Convert(Bool_t toGMT=kFALSE) const
Convert fDatime from TDatime format to the standard time_t format.
Definition TDatime.cxx:182
const char * AsString() const
Return the date & time as a string (ctime() format).
Definition TDatime.cxx:102
virtual Int_t GetValue(const char *name, Int_t dflt) const
Returns the integer value for a resource.
Definition TEnv.cxx:491
A TGraph is an object made of two arrays X and Y with npoints each.
Definition TGraph.h:41
virtual TCanvasImp * CreateCanvasImp(TCanvas *c, const char *title, UInt_t width, UInt_t height)
Create a batch version of TCanvasImp.
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:58
Option_t * GetOption() const
void Reset()
A doubly linked list.
Definition TList.h:38
virtual void Add(TObject *obj)
Definition TList.h:81
virtual TObject * Remove(TObject *obj)
Remove object from the list.
Definition TList.cxx:822
virtual TObject * FindObject(const char *name) const
Find an object in this list using its name.
Definition TList.cxx:578
virtual const char * GetName() const
Returns name of object.
Definition TNamed.h:47
An array of TObjects.
Definition TObjArray.h:31
void Add(TObject *obj)
Definition TObjArray.h:68
Int_t GetEntries() const
Return the number of objects in array (i.e.
Mother of all ROOT objects.
Definition TObject.h:41
virtual void Clear(Option_t *="")
Definition TObject.h:119
virtual const char * GetName() const
Returns name of object.
Definition TObject.cxx:429
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition TObject.h:201
virtual TObject * Clone(const char *newname="") const
Make a clone of an object using the Streamer facility.
Definition TObject.cxx:216
virtual void Execute(const char *method, const char *params, Int_t *error=0)
Execute method on this object with the given parameter string, e.g.
Definition TObject.cxx:349
R__ALWAYS_INLINE Bool_t IsOnHeap() const
Definition TObject.h:152
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:200
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:949
virtual void ExecuteEvent(Int_t event, Int_t px, Int_t py)
Execute action corresponding to an event at (px,py).
Definition TObject.cxx:383
virtual char * GetObjectInfo(Int_t px, Int_t py) const
Returns string containing info about the object at position (px,py).
Definition TObject.cxx:458
virtual void Delete(Option_t *option="")
Delete this object.
Definition TObject.cxx:241
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:766
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition TObject.cxx:515
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:963
virtual const char * GetTitle() const
Returns title of object.
Definition TObject.cxx:473
virtual void Pop()
Pop on object drawn in a pad to the top of the display list.
Definition TObject.cxx:600
void ResetBit(UInt_t f)
Definition TObject.h:200
@ kNoContextMenu
if object does not want context menu
Definition TObject.h:69
@ kMustCleanup
if object destructor must call RecursiveRemove()
Definition TObject.h:64
virtual void Info(const char *method, const char *msgfmt,...) const
Issue info message.
Definition TObject.cxx:937
Implement TVirtualPadPainter which abstracts painting operations.
Definition TPadPainter.h:26
The most important graphics class in the ROOT system.
Definition TPad.h:26
void SetGridx(Int_t value=1) override
Definition TPad.h:327
Short_t GetBorderMode() const override
Definition TPad.h:193
void SetBorderSize(Short_t bordersize) override
Definition TPad.h:317
Int_t GetTicky() const override
Definition TPad.h:233
void PaintBorder(Color_t color, Bool_t tops)
Paint the pad border.
Definition TPad.cxx:3515
void ResizePad(Option_t *option="") override
Compute pad conversion coefficients.
Definition TPad.cxx:5490
void SetGrid(Int_t valuex=1, Int_t valuey=1) override
Definition TPad.h:326
void SetTickx(Int_t value=1) override
Definition TPad.h:347
void SavePrimitive(std::ostream &out, Option_t *option="") override
Save primitives in this pad on the C++ source file out.
Definition TPad.cxx:5715
Bool_t GetGridx() const override
Definition TPad.h:229
Double_t fX2
X of upper X coordinate.
Definition TPad.h:34
void SetLogz(Int_t value=1) override
Set Lin/Log scale for Z.
Definition TPad.cxx:5978
Double_t GetY2() const override
Definition TPad.h:237
void Close(Option_t *option="") override
Delete all primitives in pad and pad itself.
Definition TPad.cxx:999
void PaintModified() override
Traverse pad hierarchy and (re)paint only modified pads.
Definition TPad.cxx:3670
const char * GetTitle() const override
Returns title of object.
Definition TPad.h:255
Double_t fX1
X of lower X coordinate.
Definition TPad.h:32
TList * GetListOfPrimitives() const override
Definition TPad.h:239
void SetFillStyle(Style_t fstyle) override
Override TAttFill::FillStyle for TPad because we want to handle style=0 as style 4000.
Definition TPad.cxx:5941
void SetPad(const char *name, const char *title, Double_t xlow, Double_t ylow, Double_t xup, Double_t yup, Color_t color=35, Short_t bordersize=5, Short_t bordermode=-1) override
Set all pad parameters.
Definition TPad.cxx:6038
void UseCurrentStyle() override
Force a copy of current style for all objects in pad.
Definition TPad.cxx:6792
void Range(Double_t x1, Double_t y1, Double_t x2, Double_t y2) override
Set world coordinate system for the pad.
Definition TPad.cxx:5209
Double_t fY1
Y of lower Y coordinate.
Definition TPad.h:33
Int_t fGLDevice
! OpenGL off-screen pixmap identifier
Definition TPad.h:81
void Update() override
Update pad.
Definition TPad.cxx:2858
void Clear(Option_t *option="") override
Delete all pad primitives.
Definition TPad.cxx:634
Int_t GetTickx() const override
Definition TPad.h:232
Double_t fAspectRatio
ratio of w/h in case of fixed ratio
Definition TPad.h:78
void Modified(Bool_t flag=1) override
Definition TPad.h:413
void SetLogy(Int_t value=1) override
Set Lin/Log scale for Y.
Definition TPad.cxx:5967
void HighLight(Color_t col=kRed, Bool_t set=kTRUE) override
Highlight pad.
Definition TPad.cxx:2966
TCanvas * fCanvas
! Pointer to mother canvas
Definition TPad.h:102
Bool_t fFixedAspectRatio
True if fixed aspect ratio.
Definition TPad.h:100
void ls(Option_t *option="") const override
List all primitives in pad.
Definition TPad.cxx:3001
TString fTitle
Pad title.
Definition TPad.h:106
void CopyPixmaps() override
Copy the sub-pixmaps of the pad to the canvas.
Definition TPad.cxx:1071
void CopyPixmap() override
Copy the pixmap of the pad to the canvas.
Definition TPad.cxx:1057
Double_t GetY1() const override
Definition TPad.h:236
Bool_t GetGridy() const override
Definition TPad.h:230
void ExecuteEvent(Int_t event, Int_t px, Int_t py) override
Execute action corresponding to one event.
Definition TPad.cxx:1717
Int_t GetLogz() const override
Definition TPad.h:252
Short_t GetBorderSize() const override
Definition TPad.h:194
TList * fPrimitives
->List of primitives (subpads)
Definition TPad.h:103
TCanvas * GetCanvas() const override
Definition TPad.h:256
Short_t fBorderSize
pad bordersize in pixels
Definition TPad.h:93
void SetGridy(Int_t value=1) override
Definition TPad.h:328
void Paint(Option_t *option="") override
Paint all primitives in pad.
Definition TPad.cxx:3453
TString fName
Pad name.
Definition TPad.h:105
Int_t fPixmapID
! Off-screen pixmap identifier
Definition TPad.h:80
TObject * FindObject(const char *name) const override
Search if object named name is inside this pad or in pads inside this pad.
Definition TPad.cxx:2628
Color_t GetHighLightColor() const override
Get highlight color.
Definition TPad.cxx:2716
TVirtualPad * cd(Int_t subpadnumber=0) override
Set Current pad.
Definition TPad.cxx:604
Int_t GetLogy() const override
Definition TPad.h:251
void SetBorderMode(Short_t bordermode) override
Definition TPad.h:316
void SetTicks(Int_t valuex=1, Int_t valuey=1) override
Definition TPad.h:346
void SetTicky(Int_t value=1) override
Definition TPad.h:348
Double_t fY2
Y of upper Y coordinate.
Definition TPad.h:35
TVirtualPad * GetSelectedPad() const override
Get selected pad.
Definition TPad.cxx:2741
Short_t fBorderMode
Bordermode (-1=down, 0 = no border, 1=up)
Definition TPad.h:94
void SetLogx(Int_t value=1) override
Set Lin/Log scale for X.
Definition TPad.cxx:5953
Int_t GetLogx() const override
Definition TPad.h:250
Double_t GetX2() const override
Definition TPad.h:235
Double_t GetX1() const override
Definition TPad.h:234
TPad * fMother
! pointer to mother of the list
Definition TPad.h:101
const char * GetName() const override
Returns name of object.
Definition TPad.h:254
void Emit(const char *signal, const T &arg)
Activate signal with single parameter.
Definition TQObject.h:164
Bool_t Connect(const char *signal, const char *receiver_class, void *receiver, const char *slot)
Non-static method is used to connect from the signal of this object to the receiver slot.
Definition TQObject.cxx:869
static Int_t IncreaseDirLevel()
Increase the indentation level for ls().
Definition TROOT.cxx:2811
static void IndentLevel()
Functions used by ls() to indent an object hierarchy.
Definition TROOT.cxx:2819
static Int_t DecreaseDirLevel()
Decrease the indentation level for ls().
Definition TROOT.cxx:2704
Basic string class.
Definition TString.h:136
Ssiz_t Length() const
Definition TString.h:410
void ToLower()
Change string to lower-case.
Definition TString.cxx:1150
Ssiz_t Last(char c) const
Find last occurrence of a character c.
Definition TString.cxx:916
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:624
Ssiz_t Index(const char *pat, Ssiz_t i=0, ECaseCompare cmp=kExact) const
Definition TString.h:639
Double_t GetTimeOffset() const
Definition TStyle.h:260
Int_t GetOptLogy() const
Definition TStyle.h:239
Int_t GetOptStat() const
Definition TStyle.h:236
Int_t GetOptTitle() const
Definition TStyle.h:237
void SetCanvasBorderSize(Width_t size=1)
Definition TStyle.h:328
Float_t GetScreenFactor() const
Definition TStyle.h:247
Int_t GetPadTickX() const
Definition TStyle.h:208
Bool_t IsReading() const
Definition TStyle.h:282
void SetCanvasColor(Color_t color=19)
Definition TStyle.h:327
Float_t GetPadRightMargin() const
Definition TStyle.h:205
void SetCanvasBorderMode(Int_t mode=1)
Definition TStyle.h:329
void SetPalette(Int_t ncolors=kBird, Int_t *colors=0, Float_t alpha=1.)
See TColor::SetPalette.
Definition TStyle.cxx:1782
Int_t GetCanvasDefH() const
Definition TStyle.h:183
Int_t GetCanvasDefX() const
Definition TStyle.h:185
Bool_t GetPadGridY() const
Definition TStyle.h:207
Float_t GetPadLeftMargin() const
Definition TStyle.h:204
Bool_t GetCanvasPreferGL() const
Definition TStyle.h:179
Int_t GetCanvasDefY() const
Definition TStyle.h:186
Bool_t GetPadGridX() const
Definition TStyle.h:206
Int_t GetPadTickY() const
Definition TStyle.h:209
Color_t GetCanvasColor() const
Definition TStyle.h:180
Float_t GetPadBottomMargin() const
Definition TStyle.h:202
Int_t GetCanvasDefW() const
Definition TStyle.h:184
Int_t GetOptLogx() const
Definition TStyle.h:238
Int_t GetCanvasBorderMode() const
Definition TStyle.h:182
Width_t GetCanvasBorderSize() const
Definition TStyle.h:181
Int_t GetOptFit() const
Definition TStyle.h:235
Int_t GetOptLogz() const
Definition TStyle.h:240
Float_t GetPadTopMargin() const
Definition TStyle.h:203
TVirtualPS is an abstract interface to Postscript, PDF, SVG.
Definition TVirtualPS.h:30
static TVirtualPadEditor * GetPadEditor(Bool_t load=kTRUE)
Returns the pad editor dialog. Static method.
To make it possible to use GL for 2D graphic in a TPad/TCanvas.
virtual void LockPainter()
Empty definition.
static TVirtualPadPainter * PadPainter(Option_t *opt="")
Create a pad painter of specified type.
virtual void SelectDrawable(Int_t device)=0
virtual void InitPainter()
Empty definition.
TVirtualPad is an abstract base class for the Pad and Canvas classes.
Definition TVirtualPad.h:51
virtual void SetCursor(ECursor cursor)=0
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
R__ALWAYS_INLINE bool HasBeenDeleted(const TObject *obj)
Check if the TObject's memory has been deleted.
Definition TObject.h:404
Int_t Nint(T x)
Round to nearest integer. Rounds half integers to the nearest even integer.
Definition TMath.h:663
Definition first.py:1
th1 Draw()