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
67{
68 const char *defcanvas = gROOT->GetDefCanvasName();
69 TString cdef = defcanvas;
70
71 auto lc = (TList*)gROOT->GetListOfCanvases();
72 Int_t n = lc->GetSize() + 1;
73
74 while(lc->FindObject(cdef.Data()))
75 cdef.Form("%s_n%d", defcanvas, n++);
76
77 return cdef;
78}
79
80
81/** \class TCanvas
82\ingroup gpad
83
84The Canvas class.
85
86A Canvas is an area mapped to a window directly under the control of the display
87manager. A ROOT session may have several canvases open at any given time.
88
89A Canvas may be subdivided into independent graphical areas: the __Pads__.
90A canvas has a default pad which has the name of the canvas itself.
91An example of a Canvas layout is sketched in the picture below.
92
93\image html gpad_canvas.png
94
95This canvas contains two pads named P1 and P2. Both Canvas, P1 and P2 can be
96moved, grown, shrunk using the normal rules of the Display manager.
97
98Once objects have been drawn in a canvas, they can be edited/moved by pointing
99directly to them. The cursor shape is changed to suggest the type of action that
100one can do on this object. Clicking with the right mouse button on an object
101pops-up a contextmenu with a complete list of actions possible on this object.
102
103A graphical editor may be started from the canvas "View" menu under the menu
104entry "Toolbar".
105
106An interactive HELP is available by clicking on the HELP button at the top right
107of the canvas. It gives a short explanation about the canvas' menus.
108
109A canvas may be automatically divided into pads via `TPad::Divide`.
110
111At creation time, no matter if in interactive or batch mode, the constructor
112defines the size of the canvas window (including the size of the window
113manager's decoration). To define precisely the graphics area size of a canvas in
114the interactive mode, the following four lines of code should be used:
115~~~ {.cpp}
116 {
117 Double_t w = 600;
118 Double_t h = 600;
119 auto c = new TCanvas("c", "c", w, h);
120 c->SetWindowSize(w + (w - c->GetWw()), h + (h - c->GetWh()));
121 }
122~~~
123and in the batch mode simply do:
124~~~ {.cpp}
125 c->SetCanvasSize(w,h);
126~~~
127
128If the canvas size exceeds the window size, scroll bars will be added to the canvas
129This allows to display very large canvases (even bigger than the screen size). The
130Following example shows how to proceed.
131~~~ {.cpp}
132 {
133 auto c = new TCanvas("c","c");
134 c->SetCanvasSize(1500, 1500);
135 c->SetWindowSize(500, 500);
136 }
137~~~
138*/
139
140////////////////////////////////////////////////////////////////////////////////
141/// Canvas default constructor.
142
143TCanvas::TCanvas(Bool_t build) : TPad(), fDoubleBuffer(0)
144{
145 fPainter = nullptr;
146 fWindowTopX = 0;
147 fWindowTopY = 0;
148 fWindowWidth = 0;
149 fWindowHeight = 0;
150 fCw = 0;
151 fCh = 0;
152 fXsizeUser = 0;
153 fYsizeUser = 0;
156 fHighLightColor = gEnv->GetValue("Canvas.HighLightColor", kRed);
157 fEvent = -1;
158 fEventX = -1;
159 fEventY = -1;
160 fSelectedX = 0;
161 fSelectedY = 0;
163 fDrawn = kFALSE;
164 fSelected = nullptr;
165 fClickSelected = nullptr;
166 fSelectedPad = nullptr;
167 fClickSelectedPad = nullptr;
168 fPadSave = nullptr;
169 fCanvasImp = nullptr;
170 fContextMenu = nullptr;
171
173
174 if (!build || TClass::IsCallingNew() != TClass::kRealNew) {
175 Constructor();
176 } else {
177 TString cdef = GetNewCanvasName();
178
179 Constructor(cdef.Data(), cdef.Data(), 1);
180 }
181}
182
183////////////////////////////////////////////////////////////////////////////////
184/// Canvas default constructor
185
187{
188 if (gThreadXAR) {
189 void *arr[2];
190 arr[1] = this;
191 if ((*gThreadXAR)("CANV", 2, arr, nullptr)) return;
192 }
193
194 fCanvas = nullptr;
195 fCanvasID = -1;
196 fCanvasImp = nullptr;
197 fBatch = kTRUE;
199
200 fContextMenu = nullptr;
201 fSelected = nullptr;
202 fClickSelected = nullptr;
203 fSelectedPad = nullptr;
204 fClickSelectedPad = nullptr;
205 fPadSave = nullptr;
209}
210
211////////////////////////////////////////////////////////////////////////////////
212/// Create an embedded canvas, i.e. a canvas that is in a TGCanvas widget
213/// which is placed in a TGFrame. This ctor is only called via the
214/// TRootEmbeddedCanvas class.
215///
216/// If "name" starts with "gl" the canvas is ready to receive GL output.
217
218TCanvas::TCanvas(const char *name, Int_t ww, Int_t wh, Int_t winid) : TPad(), fDoubleBuffer(0)
219{
220 fCanvasImp = nullptr;
221 fPainter = nullptr;
222 Init();
223
224 fCanvasID = winid;
225 fWindowTopX = 0;
226 fWindowTopY = 0;
227 fWindowWidth = ww;
228 fWindowHeight = wh;
229 fCw = ww + 4;
230 fCh = wh +28;
231 fBatch = kFALSE;
233
234 //This is a very special ctor. A window exists already!
235 //Can create painter now.
237
238 if (fUseGL) {
239 fGLDevice = gGLManager->CreateGLContext(winid);
240 if (fGLDevice == -1)
241 fUseGL = kFALSE;
242 }
243
245 if (!fCanvasImp) return;
246
248 SetName(name);
249 Build();
250}
251
252////////////////////////////////////////////////////////////////////////////////
253/// Create a new canvas with a predefined size form.
254/// If form < 0 the menubar is not shown.
255///
256/// - form = 1 700x500 at 10,10 (set by TStyle::SetCanvasDefH,W,X,Y)
257/// - form = 2 500x500 at 20,20
258/// - form = 3 500x500 at 30,30
259/// - form = 4 500x500 at 40,40
260/// - form = 5 500x500 at 50,50
261///
262/// If "name" starts with "gl" the canvas is ready to receive GL output.
263
264TCanvas::TCanvas(const char *name, const char *title, Int_t form) : TPad(), fDoubleBuffer(0)
265{
266 fPainter = nullptr;
268
269 Constructor(name, title, form);
270}
271
272////////////////////////////////////////////////////////////////////////////////
273/// Create a new canvas with a predefined size form.
274/// If form < 0 the menubar is not shown.
275///
276/// - form = 1 700x500 at 10,10 (set by TStyle::SetCanvasDefH,W,X,Y)
277/// - form = 2 500x500 at 20,20
278/// - form = 3 500x500 at 30,30
279/// - form = 4 500x500 at 40,40
280/// - form = 5 500x500 at 50,50
281
282void TCanvas::Constructor(const char *name, const char *title, Int_t form)
283{
284 if (gThreadXAR) {
285 void *arr[6];
286 static Int_t ww = 500;
287 static Int_t wh = 500;
288 arr[1] = this; arr[2] = (void*)name; arr[3] = (void*)title; arr[4] =&ww; arr[5] = &wh;
289 if ((*gThreadXAR)("CANV", 6, arr, nullptr)) return;
290 }
291
292 Init();
293 SetBit(kMenuBar,1);
294 if (form < 0) {
295 form = -form;
296 SetBit(kMenuBar,0);
297 }
298
299 fCanvas = this;
300
301 fCanvasID = -1;
302 TCanvas *old = (TCanvas*)gROOT->GetListOfCanvases()->FindObject(name);
303 if (old && old->IsOnHeap()) {
304 Warning("Constructor","Deleting canvas with same name: %s",name);
305 delete old;
306 }
307 if (gROOT->IsBatch()) { //We are in Batch mode
309 if (form == 1) {
312 } else {
313 fWindowWidth = 500;
314 fWindowHeight = 500;
315 }
319 if (!fCanvasImp) return;
320 fBatch = kTRUE;
321 } else { //normal mode with a screen window
323 if (form < 1 || form > 5) form = 1;
324 auto factory = gROOT->IsWebDisplay() ? gBatchGuiFactory : gGuiFactory;
325 if (form == 1) {
326 UInt_t uh = UInt_t(cx*gStyle->GetCanvasDefH());
327 UInt_t uw = UInt_t(cx*gStyle->GetCanvasDefW());
328 Int_t ux = Int_t(cx*gStyle->GetCanvasDefX());
329 Int_t uy = Int_t(cx*gStyle->GetCanvasDefY());
330 fCanvasImp = factory->CreateCanvasImp(this, name, ux, uy, uw, uh);
331 }
332 fCw = 500;
333 fCh = 500;
334 if (form == 2) fCanvasImp = factory->CreateCanvasImp(this, name, 20, 20, UInt_t(cx*500), UInt_t(cx*500));
335 if (form == 3) fCanvasImp = factory->CreateCanvasImp(this, name, 30, 30, UInt_t(cx*500), UInt_t(cx*500));
336 if (form == 4) fCanvasImp = factory->CreateCanvasImp(this, name, 40, 40, UInt_t(cx*500), UInt_t(cx*500));
337 if (form == 5) fCanvasImp = factory->CreateCanvasImp(this, name, 50, 50, UInt_t(cx*500), UInt_t(cx*500));
338 if (!fCanvasImp) return;
339
340 if (!gROOT->IsBatch() && fCanvasID == -1)
342
344 fBatch = kFALSE;
345 }
346
348
349 SetName(name);
350 SetTitle(title); // requires fCanvasImp set
351 Build();
352
353 // Popup canvas
354 fCanvasImp->Show();
355}
356
357////////////////////////////////////////////////////////////////////////////////
358/// Create a new canvas at a random position.
359///
360/// \param[in] name canvas name
361/// \param[in] title canvas title
362/// \param[in] ww is the window size in pixels along X
363/// (if ww < 0 the menubar is not shown)
364/// \param[in] wh is the window size in pixels along Y
365///
366/// If "name" starts with "gl" the canvas is ready to receive GL output.
367
368TCanvas::TCanvas(const char *name, const char *title, Int_t ww, Int_t wh) : TPad(), fDoubleBuffer(0)
369{
370 fPainter = nullptr;
372
373 Constructor(name, title, ww, wh);
374}
375
376////////////////////////////////////////////////////////////////////////////////
377/// Create a new canvas at a random position.
378///
379/// \param[in] name canvas name
380/// \param[in] title canvas title
381/// \param[in] ww is the window size in pixels along X
382/// (if ww < 0 the menubar is not shown)
383/// \param[in] wh is the window size in pixels along Y
384
385void TCanvas::Constructor(const char *name, const char *title, Int_t ww, Int_t wh)
386{
387 if (gThreadXAR) {
388 void *arr[6];
389 arr[1] = this; arr[2] = (void*)name; arr[3] = (void*)title; arr[4] =&ww; arr[5] = &wh;
390 if ((*gThreadXAR)("CANV", 6, arr, nullptr)) return;
391 }
392
393 Init();
394 SetBit(kMenuBar,1);
395 if (ww < 0) {
396 ww = -ww;
397 SetBit(kMenuBar,0);
398 }
399 if (wh <= 0) {
400 Error("Constructor", "Invalid canvas height: %d",wh);
401 return;
402 }
403 fCw = ww;
404 fCh = wh;
405 fCanvasID = -1;
406 TCanvas *old = (TCanvas*)gROOT->GetListOfCanvases()->FindObject(name);
407 if (old && old->IsOnHeap()) {
408 Warning("Constructor","Deleting canvas with same name: %s",name);
409 delete old;
410 }
411 if (gROOT->IsBatch()) { //We are in Batch mode
413 fWindowWidth = ww;
414 fWindowHeight = wh;
415 fCw = ww;
416 fCh = wh;
418 if (!fCanvasImp) return;
419 fBatch = kTRUE;
420 } else {
422 auto factory = gROOT->IsWebDisplay() ? gBatchGuiFactory : gGuiFactory;
423 fCanvasImp = factory->CreateCanvasImp(this, name, UInt_t(cx*ww), UInt_t(cx*wh));
424 if (!fCanvasImp) return;
425
426 if (!gROOT->IsBatch() && fCanvasID == -1)
428
430 fBatch = kFALSE;
431 }
432
434
435 SetName(name);
436 SetTitle(title); // requires fCanvasImp set
437 Build();
438
439 // Popup canvas
440 fCanvasImp->Show();
441}
442
443////////////////////////////////////////////////////////////////////////////////
444/// Create a new canvas.
445///
446/// \param[in] name canvas name
447/// \param[in] title canvas title
448/// \param[in] wtopx,wtopy are the pixel coordinates of the top left corner of
449/// the canvas (if wtopx < 0) the menubar is not shown)
450/// \param[in] ww is the window size in pixels along X
451/// \param[in] wh is the window size in pixels along Y
452///
453/// If "name" starts with "gl" the canvas is ready to receive GL output.
454
455TCanvas::TCanvas(const char *name, const char *title, Int_t wtopx, Int_t wtopy, Int_t ww, Int_t wh)
456 : TPad(), fDoubleBuffer(0)
457{
458 fPainter = nullptr;
460
461 Constructor(name, title, wtopx, wtopy, ww, wh);
462}
463
464////////////////////////////////////////////////////////////////////////////////
465/// Create a new canvas.
466///
467/// \param[in] name canvas name
468/// \param[in] title canvas title
469/// \param[in] wtopx,wtopy are the pixel coordinates of the top left corner of
470/// the canvas (if wtopx < 0) the menubar is not shown)
471/// \param[in] ww is the window size in pixels along X
472/// \param[in] wh is the window size in pixels along Y
473
474void TCanvas::Constructor(const char *name, const char *title, Int_t wtopx,
475 Int_t wtopy, Int_t ww, Int_t wh)
476{
477 if (gThreadXAR) {
478 void *arr[8];
479 arr[1] = this; arr[2] = (void*)name; arr[3] = (void*)title;
480 arr[4] = &wtopx; arr[5] = &wtopy; arr[6] = &ww; arr[7] = &wh;
481 if ((*gThreadXAR)("CANV", 8, arr, nullptr)) return;
482 }
483
484 Init();
485 SetBit(kMenuBar,1);
486 if (wtopx < 0) {
487 wtopx = -wtopx;
488 SetBit(kMenuBar,0);
489 }
490 fCw = ww;
491 fCh = wh;
492 fCanvasID = -1;
493 TCanvas *old = (TCanvas*)gROOT->GetListOfCanvases()->FindObject(name);
494 if (old && old->IsOnHeap()) {
495 Warning("Constructor","Deleting canvas with same name: %s",name);
496 delete old;
497 }
498 if (gROOT->IsBatch()) { //We are in Batch mode
500 fWindowWidth = ww;
501 fWindowHeight = wh;
502 fCw = ww;
503 fCh = wh;
505 if (!fCanvasImp) return;
506 fBatch = kTRUE;
507 } else { //normal mode with a screen window
509 auto factory = gROOT->IsWebDisplay() ? gBatchGuiFactory : gGuiFactory;
510 fCanvasImp = factory->CreateCanvasImp(this, name, Int_t(cx*wtopx), Int_t(cx*wtopy), UInt_t(cx*ww), UInt_t(cx*wh));
511 if (!fCanvasImp) return;
512
513 if (!gROOT->IsBatch() && fCanvasID == -1)
515
517 fBatch = kFALSE;
518 }
519
521
522 SetName(name);
523 SetTitle(title); // requires fCanvasImp set
524 Build();
525
526 // Popup canvas
527 fCanvasImp->Show();
528}
529
530////////////////////////////////////////////////////////////////////////////////
531/// Initialize the TCanvas members. Called by all constructors.
532
534{
535 // Make sure the application environment exists. It is need for graphics
536 // (colors are initialized in the TApplication ctor).
537 if (!gApplication)
539
540 // Load and initialize graphics libraries if
541 // TApplication::NeedGraphicsLibs() has been called by a
542 // library static initializer.
543 if (gApplication)
544 gApplication->InitializeGraphics(gROOT->IsWebDisplay());
545
546 // Get some default from .rootrc. Used in fCanvasImp->InitWindow().
547 fHighLightColor = gEnv->GetValue("Canvas.HighLightColor", kRed);
548 SetBit(kMoveOpaque, gEnv->GetValue("Canvas.MoveOpaque", 0));
549 SetBit(kResizeOpaque, gEnv->GetValue("Canvas.ResizeOpaque", 0));
550 if (gEnv->GetValue("Canvas.ShowEventStatus", kFALSE)) SetBit(kShowEventStatus);
551 if (gEnv->GetValue("Canvas.ShowToolTips", kFALSE)) SetBit(kShowToolTips);
552 if (gEnv->GetValue("Canvas.ShowToolBar", kFALSE)) SetBit(kShowToolBar);
553 if (gEnv->GetValue("Canvas.ShowEditor", kFALSE)) SetBit(kShowEditor);
554 if (gEnv->GetValue("Canvas.AutoExec", kTRUE)) SetBit(kAutoExec);
555
556 // Fill canvas ROOT data structure
557 fXsizeUser = 0;
558 fYsizeUser = 0;
561
562 fDISPLAY = "$DISPLAY";
565 fSelected = nullptr;
566 fClickSelected = nullptr;
567 fSelectedX = 0;
568 fSelectedY = 0;
569 fSelectedPad = nullptr;
570 fClickSelectedPad= nullptr;
571 fPadSave = nullptr;
572 fEvent = -1;
573 fEventX = -1;
574 fEventY = -1;
575 fContextMenu = nullptr;
576 fDrawn = kFALSE;
577}
578
579////////////////////////////////////////////////////////////////////////////////
580/// Build a canvas. Called by all constructors.
581
583{
584 // Get window identifier
585 if (fCanvasID == -1 && fCanvasImp)
587 if (fCanvasID == -1) return;
588
589 if (fCw !=0 && fCh !=0) {
592 }
593
594 // Set Pad parameters
595 gPad = this;
596 fCanvas = this;
597 fMother = (TPad*)gPad;
598
599 if (IsBatch()) {
600 // Make sure that batch interactive canvas sizes are the same
601 fCw -= 4;
602 fCh -= 28;
603 } else if (IsWeb()) {
604 // mark canvas as batch - avoid virtualx in many places
606 } else {
607 //normal mode with a screen window
608 // Set default physical canvas attributes
609 //Should be done via gVirtualX, not via fPainter (at least now). No changes here.
610 gVirtualX->SelectWindow(fCanvasID);
611 gVirtualX->SetFillColor(1); //Set color index for fill area
612 gVirtualX->SetLineColor(1); //Set color index for lines
613 gVirtualX->SetMarkerColor(1); //Set color index for markers
614 gVirtualX->SetTextColor(1); //Set color index for text
615 // Clear workstation
616 gVirtualX->ClearWindow();
617
618 // Set Double Buffer on by default
620
621 // Get effective window parameters (with borders and menubar)
624
625 // Get effective canvas parameters without borders
626 Int_t dum1, dum2;
627 gVirtualX->GetGeometry(fCanvasID, dum1, dum2, fCw, fCh);
628
629 fContextMenu = new TContextMenu("ContextMenu");
630 }
631
632 gROOT->GetListOfCanvases()->Add(this);
633
634 if (!fPrimitives) {
635 fPrimitives = new TList;
637 SetFillStyle(1001);
649 fBorderMode=gStyle->GetCanvasBorderMode(); // do not call SetBorderMode (function redefined in TCanvas)
650 SetPad(0, 0, 1, 1);
651 Range(0, 0, 1, 1); //pad range is set by default to [0,1] in x and y
652
654 if (vpp) vpp->SelectDrawable(fPixmapID);//gVirtualX->SelectPixmap(fPixmapID); //pixmap must be selected
655 PaintBorder(GetFillColor(), kTRUE); //paint background
656 }
657
658 // transient canvases have typically no menubar and should not get
659 // by default the event status bar (if set by default)
660 if (TestBit(kMenuBar) && fCanvasImp) {
662 // ... and toolbar + editor
666 }
667}
668
669////////////////////////////////////////////////////////////////////////////////
670/// Canvas destructor
671
673{
674 Destructor();
675}
676
677////////////////////////////////////////////////////////////////////////////////
678/// Browse.
679
681{
682 Draw();
683 cd();
685}
686
687////////////////////////////////////////////////////////////////////////////////
688/// Actual canvas destructor.
689
691{
692 if (gThreadXAR) {
693 void *arr[2];
694 arr[1] = this;
695 if ((*gThreadXAR)("CDEL", 2, arr, nullptr)) return;
696 }
697
698 if (ROOT::Detail::HasBeenDeleted(this)) return;
699
701 if (!gPad) return;
702
703 Close();
704
705 //If not yet (batch mode?).
707}
708
709////////////////////////////////////////////////////////////////////////////////
710/// Set current canvas & pad. Returns the new current pad,
711/// or 0 in case of failure.
712/// See TPad::cd() for an explanation of the parameter.
713
715{
716 if (fCanvasID == -1) return nullptr;
717
718 TPad::cd(subpadnumber);
719
720 // in case doublebuffer is off, draw directly onto display window
721 if (!IsBatch()) {
722 if (!fDoubleBuffer)
723 gVirtualX->SelectWindow(fCanvasID);//Ok, does not matter for glpad.
724 }
725 return gPad;
726}
727
728////////////////////////////////////////////////////////////////////////////////
729/// Remove all primitives from the canvas.
730/// If option "D" is specified, direct sub-pads are cleared but not deleted.
731/// This option is not recursive, i.e. pads in direct sub-pads are deleted.
732
734{
735 if (fCanvasID == -1) return;
736
738
739 TString opt = option;
740 opt.ToLower();
741 if (opt.Contains("d")) {
742 // clear subpads, but do not delete pads in case the canvas
743 // has been divided (note: option "D" is propagated so could cause
744 // conflicts for primitives using option "D" for something else)
745 if (fPrimitives) {
746 TIter next(fPrimitives);
747 TObject *obj;
748 while ((obj=next())) {
749 obj->Clear(option);
750 }
751 }
752 } else {
753 //default, clear everything in the canvas. Subpads are deleted
754 TPad::Clear(option); //Remove primitives from pad
755 }
756
757 fSelected = nullptr;
758 fClickSelected = nullptr;
759 fSelectedPad = nullptr;
760 fClickSelectedPad = nullptr;
761}
762
763////////////////////////////////////////////////////////////////////////////////
764/// Emit pad Cleared signal.
765
767{
768 Emit("Cleared(TVirtualPad*)", (Longptr_t)pad);
769}
770
771////////////////////////////////////////////////////////////////////////////////
772/// Emit Closed signal.
773
775{
776 Emit("Closed()");
777}
778
779////////////////////////////////////////////////////////////////////////////////
780/// Close canvas.
781///
782/// Delete window/pads data structure
783
785{
786 auto padsave = gPad;
787 TCanvas *cansave = padsave ? padsave->GetCanvas() : nullptr;
788
789 if (fCanvasID != -1) {
790
791 if ((!gROOT->IsLineProcessing()) && (!gVirtualX->IsCmdThread())) {
792 gInterpreter->Execute(this, IsA(), "Close", option);
793 return;
794 }
795
797
799
800 cd();
802
803 if (!IsBatch()) {
804 gVirtualX->SelectWindow(fCanvasID); //select current canvas
805
807
808 if (fCanvasImp)
809 fCanvasImp->Close();
810 }
811 fCanvasID = -1;
812 fBatch = kTRUE;
813
814 gROOT->GetListOfCanvases()->Remove(this);
815
816 // Close actual window on screen
818 }
819
820 if (cansave == this) {
821 gPad = (TCanvas *) gROOT->GetListOfCanvases()->First();
822 } else {
823 gPad = padsave;
824 }
825
826 Closed();
827}
828
829////////////////////////////////////////////////////////////////////////////////
830/// Copy the canvas pixmap of the pad to the canvas.
831
833{
834 if (!IsBatch()) {
835 CopyPixmap();
837 }
838}
839
840////////////////////////////////////////////////////////////////////////////////
841/// Draw a canvas.
842/// If a canvas with the name is already on the screen, the canvas is repainted.
843/// This function is useful when a canvas object has been saved in a Root file.
844/// One can then do:
845/// ~~~ {.cpp}
846/// Root > TFile::Open("file.root");
847/// Root > canvas->Draw();
848/// ~~~
849
851{
852 // Load and initialize graphics libraries if
853 // TApplication::NeedGraphicsLibs() has been called by a
854 // library static initializer.
855 if (gApplication)
856 gApplication->InitializeGraphics(gROOT->IsWebDisplay());
857
858 fDrawn = kTRUE;
859
860 TCanvas *old = (TCanvas*)gROOT->GetListOfCanvases()->FindObject(GetName());
861 if (old == this) {
862 Paint();
863 return;
864 }
865 if (old) { gROOT->GetListOfCanvases()->Remove(old); delete old;}
866
867 if (fWindowWidth == 0) {
868 if (fCw !=0) fWindowWidth = fCw+4;
869 else fWindowWidth = 800;
870 }
871 if (fWindowHeight == 0) {
872 if (fCh !=0) fWindowHeight = fCh+28;
873 else fWindowHeight = 600;
874 }
875 if (gROOT->IsBatch()) { //We are in Batch mode
877 if (!fCanvasImp) return;
878 fBatch = kTRUE;
879
880 } else { //normal mode with a screen window
881 auto factory = gROOT->IsWebDisplay() ? gBatchGuiFactory : gGuiFactory;
882 fCanvasImp = factory->CreateCanvasImp(this, GetName(), fWindowTopX, fWindowTopY,
884 if (!fCanvasImp) return;
886 }
887 Build();
888 ResizePad();
890 fCanvasImp->Show();
891 Modified();
892}
893
894////////////////////////////////////////////////////////////////////////////////
895/// Draw a clone of this canvas
896/// A new canvas is created that is a clone of this canvas
897
899{
900 TCanvas *newCanvas = (TCanvas*)Clone();
901 newCanvas->SetName();
902
903 newCanvas->Draw(option);
904 newCanvas->Update();
905 return newCanvas;
906}
907
908////////////////////////////////////////////////////////////////////////////////
909/// Draw a clone of this canvas into the current pad
910/// In an interactive session, select the destination/current pad
911/// with the middle mouse button, then point to the canvas area to select
912/// the canvas context menu item DrawClonePad.
913/// Note that the original canvas may have subpads.
914
916{
917 auto padsav = gPad;
918 auto selpad = gROOT->GetSelectedPad();
919 auto pad = padsav;
920 if (pad == this) pad = selpad;
921 if (!padsav || !pad || pad == this) {
922 TCanvas *newCanvas = (TCanvas*)DrawClone();
924 return newCanvas;
925 }
926 if (fCanvasID == -1) {
927 auto factory = gROOT->IsWebDisplay() ? gBatchGuiFactory : gGuiFactory;
928 fCanvasImp = factory->CreateCanvasImp(this, GetName(), fWindowTopX, fWindowTopY,
930 if (!fCanvasImp) return nullptr;
933 }
934 this->cd();
935 //copy pad attributes
936 pad->Range(fX1,fY1,fX2,fY2);
937 pad->SetTickx(GetTickx());
938 pad->SetTicky(GetTicky());
939 pad->SetGridx(GetGridx());
940 pad->SetGridy(GetGridy());
941 pad->SetLogx(GetLogx());
942 pad->SetLogy(GetLogy());
943 pad->SetLogz(GetLogz());
944 pad->SetBorderSize(GetBorderSize());
945 pad->SetBorderMode(GetBorderMode());
946 TAttLine::Copy((TAttLine&)*pad);
947 TAttFill::Copy((TAttFill&)*pad);
948 TAttPad::Copy((TAttPad&)*pad);
949
950 //copy primitives
952 while (auto obj = next()) {
953 pad->cd();
954 auto clone = obj->Clone();
955 pad->GetListOfPrimitives()->Add(clone, next.GetOption());
956 }
957 pad->ResizePad();
958 pad->Modified();
959 pad->Update();
960 if (padsav) padsav->cd();
961 return nullptr;
962}
963
964////////////////////////////////////////////////////////////////////////////////
965/// Report name and title of primitive below the cursor.
966///
967/// This function is called when the option "Event Status"
968/// in the canvas menu "Options" is selected.
969
971{
972 const Int_t kTMAX=256;
973 static char atext[kTMAX];
974
975 if (!TestBit(kShowEventStatus) || !selected) return;
976
977 if (!fCanvasImp) return; //this may happen when closing a TAttCanvas
978
980
981 fCanvasImp->SetStatusText(selected->GetTitle(),0);
982 fCanvasImp->SetStatusText(selected->GetName(),1);
983 if (event == kKeyPress)
984 snprintf(atext, kTMAX, "%c", (char) px);
985 else
986 snprintf(atext, kTMAX, "%d,%d", px, py);
987 fCanvasImp->SetStatusText(atext,2);
988
989 // Show date/time if TimeDisplay is selected
990 TAxis *xaxis = NULL;
991 if ( selected->InheritsFrom("TH1") )
992 xaxis = ((TH1*)selected)->GetXaxis();
993 else if ( selected->InheritsFrom("TGraph") )
994 xaxis = ((TGraph*)selected)->GetXaxis();
995 else if ( selected->InheritsFrom("TAxis") )
996 xaxis = (TAxis*)selected;
997 if ( xaxis != NULL && xaxis->GetTimeDisplay()) {
998 TString objinfo = selected->GetObjectInfo(px,py);
999 // check if user has overwritten GetObjectInfo and altered
1000 // the default text from TObject::GetObjectInfo "x=.. y=.."
1001 if (objinfo.Contains("x=") && objinfo.Contains("y=") ) {
1002 UInt_t toff = 0;
1003 TString time_format(xaxis->GetTimeFormat());
1004 // TimeFormat may contain offset: %F2000-01-01 00:00:00
1005 Int_t idF = time_format.Index("%F");
1006 if (idF>=0) {
1007 Int_t lnF = time_format.Length();
1008 // minimal check for correct format
1009 if (lnF - idF == 21) {
1010 time_format = time_format(idF+2, lnF);
1011 TDatime dtoff(time_format);
1012 toff = dtoff.Convert();
1013 }
1014 } else {
1015 toff = (UInt_t)gStyle->GetTimeOffset();
1016 }
1017 TDatime dt((UInt_t)gPad->AbsPixeltoX(px) + toff);
1018 snprintf(atext, kTMAX, "%s, y=%g",
1019 dt.AsSQLString(),gPad->AbsPixeltoY(py));
1020 fCanvasImp->SetStatusText(atext,3);
1021 return;
1022 }
1023 }
1024 // default
1025 fCanvasImp->SetStatusText(selected->GetObjectInfo(px,py),3);
1026}
1027
1028////////////////////////////////////////////////////////////////////////////////
1029/// Get editor bar.
1030
1032{
1034}
1035
1036////////////////////////////////////////////////////////////////////////////////
1037/// Embedded a canvas into a TRootEmbeddedCanvas. This method is only called
1038/// via TRootEmbeddedCanvas::AdoptCanvas.
1039
1041{
1042 // If fCanvasImp already exists, no need to go further.
1043 if(fCanvasImp) return;
1044
1045 fCanvasID = winid;
1046 fWindowTopX = 0;
1047 fWindowTopY = 0;
1048 fWindowWidth = ww;
1049 fWindowHeight = wh;
1050 fCw = ww;
1051 fCh = wh;
1052 fBatch = kFALSE;
1053 fUpdating = kFALSE;
1054
1056 if (!fCanvasImp) return;
1057 Build();
1058 Resize();
1059}
1060
1061////////////////////////////////////////////////////////////////////////////////
1062/// Generate kMouseEnter and kMouseLeave events depending on the previously
1063/// selected object and the currently selected object. Does nothing if the
1064/// selected object does not change.
1065
1066void TCanvas::EnterLeave(TPad *prevSelPad, TObject *prevSelObj)
1067{
1068 if (prevSelObj == fSelected) return;
1069
1070 TContext ctxt(kFALSE);
1071 Int_t sevent = fEvent;
1072
1073 if (prevSelObj) {
1074 gPad = prevSelPad;
1075 prevSelObj->ExecuteEvent(kMouseLeave, fEventX, fEventY);
1077 RunAutoExec();
1078 ProcessedEvent(kMouseLeave, fEventX, fEventY, prevSelObj); // emit signal
1079 }
1080
1082
1083 if (fSelected) {
1086 RunAutoExec();
1088 }
1089
1090 fEvent = sevent;
1091}
1092
1093////////////////////////////////////////////////////////////////////////////////
1094/// Execute action corresponding to one event.
1095///
1096/// This member function must be implemented to realize the action
1097/// corresponding to the mouse click on the object in the canvas
1098///
1099/// Only handle mouse motion events in TCanvas, all other events are
1100/// ignored for the time being
1101
1103{
1104 if (gROOT->GetEditorMode()) {
1106 return;
1107 }
1108
1109 switch (event) {
1110
1111 case kMouseMotion:
1113 break;
1114 }
1115}
1116
1117////////////////////////////////////////////////////////////////////////////////
1118/// Turn rubberband feedback mode on or off.
1119
1121{
1122 if (set) {
1123 SetDoubleBuffer(0); // turn off double buffer mode
1124 gVirtualX->SetDrawMode(TVirtualX::kInvert); // set the drawing mode to XOR mode
1125 } else {
1126 SetDoubleBuffer(1); // turn on double buffer mode
1127 gVirtualX->SetDrawMode(TVirtualX::kCopy); // set drawing mode back to normal (copy) mode
1128 }
1129}
1130
1131////////////////////////////////////////////////////////////////////////////////
1132/// Flush canvas buffers.
1133
1135{
1136 if ((fCanvasID == -1) || IsWeb()) return;
1137
1138 TContext ctxt(this, kTRUE);
1139 if (!IsBatch()) {
1140 if (!UseGL() || fGLDevice == -1) {
1141 gVirtualX->SelectWindow(fCanvasID);
1142 gPad = ctxt.GetSaved(); //don't do cd() because than also the pixmap is changed
1143 CopyPixmaps();
1144 gVirtualX->UpdateWindow(1);
1145 } else {
1146 TVirtualPS *tvps = gVirtualPS;
1147 gVirtualPS = nullptr;
1148 gGLManager->MakeCurrent(fGLDevice);
1150 Paint();
1151 if (ctxt.GetSaved() && ctxt.GetSaved()->GetCanvas() == this) {
1152 ctxt.GetSaved()->cd();
1153 ctxt.GetSaved()->HighLight(ctxt.GetSaved()->GetHighLightColor());
1154 //cd();
1155 }
1157 gGLManager->Flush(fGLDevice);
1158 gVirtualPS = tvps;
1159 }
1160 }
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 = 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 = nullptr;
1268 fSelectedPad = nullptr;
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();
1307 Bool_t resize = kFALSE;
1309 resize = ((TBox*)fSelected)->IsBeingResized();
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
1493{
1495 std::cout <<"Canvas Name=" <<GetName()<<" Title="<<GetTitle()<<" Option="<<option<<std::endl;
1499}
1500
1501////////////////////////////////////////////////////////////////////////////////
1502/// Static function to build a default canvas.
1503
1505{
1506 TString cdef = GetNewCanvasName();
1507
1508 TCanvas *c = new TCanvas(cdef.Data(), cdef.Data(), 1);
1509
1510 ::Info("TCanvas::MakeDefCanvas"," created default TCanvas with name %s", cdef.Data());
1511 return c;
1512}
1513
1514////////////////////////////////////////////////////////////////////////////////
1515/// Set option to move objects/pads in a canvas.
1516///
1517/// - set = 1 (default) graphics objects are moved in opaque mode
1518/// - set = 0 only the outline of objects is drawn when moving them
1519///
1520/// The option opaque produces the best effect. It requires however a
1521/// a reasonably fast workstation or response time.
1522
1524{
1525 SetBit(kMoveOpaque,set);
1526}
1527
1528////////////////////////////////////////////////////////////////////////////////
1529/// Paint canvas.
1530
1532{
1533 if (fCanvas)
1535}
1536
1537////////////////////////////////////////////////////////////////////////////////
1538/// Prepare for pick, call TPad::Pick() and when selected object
1539/// is different from previous then emit Picked() signal.
1540
1541TPad *TCanvas::Pick(Int_t px, Int_t py, TObject *prevSelObj)
1542{
1543 TObjLink *pickobj = nullptr;
1544
1545 fSelected = nullptr;
1546 fSelectedOpt = "";
1547 fSelectedPad = nullptr;
1548
1549 TPad *pad = Pick(px, py, pickobj);
1550 if (!pad) return nullptr;
1551
1552 if (!pickobj) {
1553 fSelected = pad;
1554 fSelectedOpt = "";
1555 } else {
1556 if (!fSelected) { // can be set via TCanvas::SetSelected()
1557 fSelected = pickobj->GetObject();
1558 fSelectedOpt = pickobj->GetOption();
1559 }
1560 }
1561 fSelectedPad = pad;
1562
1563 if (fSelected != prevSelObj)
1564 Picked(fSelectedPad, fSelected, fEvent); // emit signal
1565
1566 if ((fEvent == kButton1Down) || (fEvent == kButton2Down) || (fEvent == kButton3Down)) {
1570 Selected(fSelectedPad, fSelected, fEvent); // emit signal
1571 fSelectedX = px;
1572 fSelectedY = py;
1573 }
1574 }
1575 return pad;
1576}
1577
1578////////////////////////////////////////////////////////////////////////////////
1579/// Emit Picked() signal.
1580
1582{
1583 Longptr_t args[3];
1584
1585 args[0] = (Longptr_t) pad;
1586 args[1] = (Longptr_t) obj;
1587 args[2] = event;
1588
1589 Emit("Picked(TPad*,TObject*,Int_t)", args);
1590}
1591
1592////////////////////////////////////////////////////////////////////////////////
1593/// Emit Highlighted() signal.
1594///
1595/// - pad is pointer to pad with highlighted histogram or graph
1596/// - obj is pointer to highlighted histogram or graph
1597/// - x is highlighted x bin for 1D histogram or highlighted x-th point for graph
1598/// - y is highlighted y bin for 2D histogram (for 1D histogram or graph not in use)
1599
1601{
1602 Longptr_t args[4];
1603
1604 args[0] = (Longptr_t) pad;
1605 args[1] = (Longptr_t) obj;
1606 args[2] = x;
1607 args[3] = y;
1608
1609 Emit("Highlighted(TVirtualPad*,TObject*,Int_t,Int_t)", args);
1610}
1611
1612////////////////////////////////////////////////////////////////////////////////
1613/// This is "simplification" for function TCanvas::Connect with Highlighted
1614/// signal for specific slot.
1615///
1616/// Slot has to be defined "UserFunction(TVirtualPad *pad, TObject *obj, Int_t x, Int_t y)"
1617/// all parameters of UserFunction are taken from TCanvas::Highlighted
1618
1619void TCanvas::HighlightConnect(const char *slot)
1620{
1621 Connect("Highlighted(TVirtualPad*,TObject*,Int_t,Int_t)", nullptr, nullptr, slot);
1622}
1623
1624////////////////////////////////////////////////////////////////////////////////
1625/// Emit Selected() signal.
1626
1628{
1629 Longptr_t args[3];
1630
1631 args[0] = (Longptr_t) pad;
1632 args[1] = (Longptr_t) obj;
1633 args[2] = event;
1634
1635 Emit("Selected(TVirtualPad*,TObject*,Int_t)", args);
1636}
1637
1638////////////////////////////////////////////////////////////////////////////////
1639/// Emit ProcessedEvent() signal.
1640
1642{
1643 Longptr_t args[4];
1644
1645 args[0] = event;
1646 args[1] = x;
1647 args[2] = y;
1648 args[3] = (Longptr_t) obj;
1649
1650 Emit("ProcessedEvent(Int_t,Int_t,Int_t,TObject*)", args);
1651}
1652
1653////////////////////////////////////////////////////////////////////////////////
1654/// Recompute canvas parameters following a X11 Resize.
1655
1657{
1658 if (fCanvasID == -1) return;
1659
1660 if ((!gROOT->IsLineProcessing()) && (!gVirtualX->IsCmdThread())) {
1661 gInterpreter->Execute(this, IsA(), "Resize", "");
1662 return;
1663 }
1664
1666
1667 TContext ctxt(this, kTRUE);
1668
1669 if (!IsBatch()) {
1670 gVirtualX->SelectWindow(fCanvasID); //select current canvas
1671 gVirtualX->ResizeWindow(fCanvasID); //resize canvas and off-screen buffer
1672
1673 // Get effective window parameters including menubar and borders
1676
1677 // Get effective canvas parameters without borders
1678 Int_t dum1, dum2;
1679 gVirtualX->GetGeometry(fCanvasID, dum1, dum2, fCw, fCh);
1680 }
1681
1682 if (fXsizeUser && fYsizeUser) {
1683 UInt_t nwh = fCh;
1684 UInt_t nww = fCw;
1686 if (rxy < 1) {
1687 UInt_t twh = UInt_t(Double_t(fCw)/rxy);
1688 if (twh > fCh)
1689 nww = UInt_t(Double_t(fCh)*rxy);
1690 else
1691 nwh = twh;
1692 if (nww > fCw) {
1693 nww = fCw; nwh = twh;
1694 }
1695 if (nwh > fCh) {
1696 nwh = fCh; nww = UInt_t(Double_t(fCh)/rxy);
1697 }
1698 } else {
1699 UInt_t twh = UInt_t(Double_t(fCw)*rxy);
1700 if (twh > fCh)
1701 nwh = UInt_t(Double_t(fCw)/rxy);
1702 else
1703 nww = twh;
1704 if (nww > fCw) {
1705 nww = fCw; nwh = twh;
1706 }
1707 if (nwh > fCh) {
1708 nwh = fCh; nww = UInt_t(Double_t(fCh)*rxy);
1709 }
1710 }
1711 fCw = nww;
1712 fCh = nwh;
1713 }
1714
1715 if (fCw < fCh) {
1718 }
1719 else {
1722 }
1723
1724//*-*- Loop on all pads to recompute conversion coefficients
1726}
1727
1728
1729////////////////////////////////////////////////////////////////////////////////
1730/// Raise canvas window
1731
1733{
1734 if (fCanvasImp)
1736}
1737
1738////////////////////////////////////////////////////////////////////////////////
1739/// Set option to resize objects/pads in a canvas.
1740///
1741/// - set = 1 (default) graphics objects are resized in opaque mode
1742/// - set = 0 only the outline of objects is drawn when resizing them
1743///
1744/// The option opaque produces the best effect. It requires however a
1745/// a reasonably fast workstation or response time.
1746
1748{
1749 SetBit(kResizeOpaque,set);
1750}
1751
1752////////////////////////////////////////////////////////////////////////////////
1753/// Execute the list of TExecs in the current pad.
1754
1756{
1757 if (!TestBit(kAutoExec))
1758 return;
1759 if (gPad)
1760 ((TPad*)gPad)->AutoExec();
1761}
1762
1763////////////////////////////////////////////////////////////////////////////////
1764/// Save primitives in this canvas in C++ macro file with GUI.
1765
1766void TCanvas::SavePrimitive(std::ostream &out, Option_t *option /*= ""*/)
1767{
1768 // Write canvas options (in $TROOT or $TStyle)
1769 if (gStyle->GetOptFit()) {
1770 out<<" gStyle->SetOptFit(1);"<<std::endl;
1771 }
1772 if (!gStyle->GetOptStat()) {
1773 out<<" gStyle->SetOptStat(0);"<<std::endl;
1774 }
1775 if (!gStyle->GetOptTitle()) {
1776 out<<" gStyle->SetOptTitle(0);"<<std::endl;
1777 }
1778 if (gROOT->GetEditHistograms()) {
1779 out<<" gROOT->SetEditHistograms();"<<std::endl;
1780 }
1781 if (GetShowEventStatus()) {
1782 out<<" "<<GetName()<<"->ToggleEventStatus();"<<std::endl;
1783 }
1784 if (GetShowToolTips()) {
1785 out<<" "<<GetName()<<"->ToggleToolTips();"<<std::endl;
1786 }
1787 if (GetShowToolBar()) {
1788 out<<" "<<GetName()<<"->ToggleToolBar();"<<std::endl;
1789 }
1790 if (GetHighLightColor() != 5) {
1792 out<<" "<<GetName()<<"->SetHighLightColor(ci);" << std::endl;
1793 else
1794 out<<" "<<GetName()<<"->SetHighLightColor("<<GetHighLightColor()<<");"<<std::endl;
1795 }
1796
1797 // Now recursively scan all pads of this canvas
1798 cd();
1800}
1801
1802////////////////////////////////////////////////////////////////////////////////
1803/// Save primitives in this canvas as a C++ macro file.
1804/// This function loops on all the canvas primitives and for each primitive
1805/// calls the object SavePrimitive function.
1806/// When outputting floating point numbers, the default precision is 7 digits.
1807/// The precision can be changed (via system.rootrc) by changing the value
1808/// of the environment variable "Canvas.SavePrecision"
1809
1811{
1812 // reset bit TClass::kClassSaved for all classes
1813 TIter next(gROOT->GetListOfClasses());
1814 TClass *cl;
1815 while((cl = (TClass*)next())) {
1817 }
1818
1819 char quote = '"';
1820 std::ofstream out;
1821 TString fname;
1822 const char *cname = GetName();
1823 Bool_t invalid = kFALSE;
1824 // if filename is given, open this file, otherwise create a file
1825 // with a name equal to the canvasname.C
1826 if (filename && (strlen(filename) > 0)) {
1827 fname = filename;
1828 } else {
1829 fname = cname;
1830 fname = fname.Strip(TString::kBoth);
1831 if (fname.IsNull()) {
1832 invalid = kTRUE;
1833 fname = "c1";
1834 }
1835 fname.Append(".C");
1836 }
1837 out.open(fname.Data(), std::ios::out);
1838 if (!out.good()) {
1839 Error("SaveSource", "Cannot open file: %s", fname.Data());
1840 return;
1841 }
1842
1843 //set precision
1844 Int_t precision = gEnv->GetValue("Canvas.SavePrecision",7);
1845 out.precision(precision);
1846
1847 // Write macro header and date/time stamp
1848 TDatime t;
1850 Int_t topx,topy;
1851 UInt_t w, h;
1852 if (!fCanvasImp) {
1853 Error("SaveSource", "Cannot open TCanvas");
1854 return;
1855 }
1856 UInt_t editorWidth = fCanvasImp->GetWindowGeometry(topx,topy,w,h);
1857 w = UInt_t((fWindowWidth - editorWidth)/cx);
1858 h = UInt_t((fWindowHeight)/cx);
1859 topx = GetWindowTopX();
1860 topy = GetWindowTopY();
1861
1862 if (w == 0) {
1863 w = GetWw()+4; h = GetWh()+4;
1864 topx = 1; topy = 1;
1865 }
1866
1867 TString mname = fname;
1868 out << R"CODE(#ifdef __CLING__
1869#pragma cling optimize(0)
1870#endif
1871)CODE";
1872 Int_t p = mname.Last('.');
1873 Int_t s = mname.Last('/')+1;
1874
1875 // A named macro is generated only if the function name is valid. If not, the
1876 // macro is unnamed.
1877 TString first(mname(s,s+1));
1878 if (!first.IsDigit()) out <<"void " << mname(s,p-s) << "()" << std::endl;
1879
1880 out <<"{"<<std::endl;
1881 out <<"//=========Macro generated from canvas: "<<GetName()<<"/"<<GetTitle()<<std::endl;
1882 out <<"//========= ("<<t.AsString()<<") by ROOT version "<<gROOT->GetVersion()<<std::endl;
1883
1885 out <<std::endl<<" gStyle->SetCanvasPreferGL(kTRUE);"<<std::endl<<std::endl;
1886
1887 // Write canvas parameters (TDialogCanvas case)
1889 out<<" "<<ClassName()<<" *"<<cname<<" = new "<<ClassName()<<"("<<quote<<GetName()
1890 <<quote<<", "<<quote<<GetTitle()<<quote<<","<<w<<","<<h<<");"<<std::endl;
1891 } else {
1892 // Write canvas parameters (TCanvas case)
1893 out<<" TCanvas *"<<cname<<" = new TCanvas("<<quote<<GetName()<<quote<<", "<<quote<<GetTitle()
1894 <<quote;
1895 if (!HasMenuBar())
1896 out<<",-"<<topx<<","<<topy<<","<<w<<","<<h<<");"<<std::endl;
1897 else
1898 out<<","<<topx<<","<<topy<<","<<w<<","<<h<<");"<<std::endl;
1899 }
1900 // Write canvas options (in $TROOT or $TStyle)
1901 if (gStyle->GetOptFit()) {
1902 out<<" gStyle->SetOptFit(1);"<<std::endl;
1903 }
1904 if (!gStyle->GetOptStat()) {
1905 out<<" gStyle->SetOptStat(0);"<<std::endl;
1906 }
1907 if (!gStyle->GetOptTitle()) {
1908 out<<" gStyle->SetOptTitle(0);"<<std::endl;
1909 }
1910 if (gROOT->GetEditHistograms()) {
1911 out<<" gROOT->SetEditHistograms();"<<std::endl;
1912 }
1913 if (GetShowEventStatus()) {
1914 out<<" "<<GetName()<<"->ToggleEventStatus();"<<std::endl;
1915 }
1916 if (GetShowToolTips()) {
1917 out<<" "<<GetName()<<"->ToggleToolTips();"<<std::endl;
1918 }
1919 if (GetHighLightColor() != 5) {
1921 out<<" "<<GetName()<<"->SetHighLightColor(ci);" << std::endl;
1922 else
1923 out<<" "<<GetName()<<"->SetHighLightColor("<<GetHighLightColor()<<");"<<std::endl;
1924 }
1925
1926 // Now recursively scan all pads of this canvas
1927 cd();
1928 if (invalid) SetName("c1");
1930 // Write canvas options related to pad editor
1931 out<<" "<<GetName()<<"->SetSelected("<<GetName()<<");"<<std::endl;
1932 if (GetShowToolBar()) {
1933 out<<" "<<GetName()<<"->ToggleToolBar();"<<std::endl;
1934 }
1935 if (invalid) SetName(" ");
1936
1937 out <<"}"<<std::endl;
1938 out.close();
1939 Info("SaveSource","C++ Macro file: %s has been generated", fname.Data());
1940
1941 // reset bit TClass::kClassSaved for all classes
1942 next.Reset();
1943 while((cl = (TClass*)next())) {
1945 }
1946}
1947
1948////////////////////////////////////////////////////////////////////////////////
1949/// Toggle batch mode. However, if the canvas is created without a window
1950/// then batch mode always stays set.
1951
1952void TCanvas::SetBatch(Bool_t batch)
1953{
1954 if (gROOT->IsBatch() || IsWeb())
1955 fBatch = kTRUE;
1956 else
1957 fBatch = batch;
1958}
1959
1960////////////////////////////////////////////////////////////////////////////////
1961/// Set Width and Height of canvas to ww and wh respectively. If ww and/or wh
1962/// are greater than the current canvas window a scroll bar is automatically
1963/// generated. Use this function to zoom in a canvas and navigate via
1964/// the scroll bars. The Width and Height in this method are different from those
1965/// given in the TCanvas constructors where these two dimension include the size
1966/// of the window decoration whereas they do not in this method.
1967
1969{
1970 if (fCanvasImp) {
1971 fCanvasImp->SetCanvasSize(ww, wh);
1972 fCw = ww;
1973 fCh = wh;
1974 TContext ctxt(this, kTRUE);
1975 ResizePad();
1976 }
1977}
1978
1979////////////////////////////////////////////////////////////////////////////////
1980/// Set cursor.
1981
1983{
1984 if (IsBatch()) return;
1985 gVirtualX->SetCursor(fCanvasID, cursor);
1986}
1987
1988////////////////////////////////////////////////////////////////////////////////
1989/// Set Double Buffer On/Off.
1990
1992{
1993 if (IsBatch()) return;
1995 gVirtualX->SetDoubleBuffer(fCanvasID, mode);
1996
1997 // depending of the buffer mode set the drawing window to either
1998 // the canvas pixmap or to the canvas on-screen window
1999 if (fDoubleBuffer) {
2001 } else
2003}
2004
2005////////////////////////////////////////////////////////////////////////////////
2006/// Fix canvas aspect ratio to current value if fixed is true.
2007
2009{
2010 if (fixed) {
2011 if (!fFixedAspectRatio) {
2012 if (fCh != 0)
2014 else {
2015 Error("SetAspectRatio", "cannot fix aspect ratio, height of canvas is 0");
2016 return;
2017 }
2019 }
2020 } else {
2022 fAspectRatio = 0;
2023 }
2024}
2025
2026////////////////////////////////////////////////////////////////////////////////
2027/// If isfolder=kTRUE, the canvas can be browsed like a folder
2028/// by default a canvas is not browsable.
2029
2030void TCanvas::SetFolder(Bool_t isfolder)
2031{
2032 fgIsFolder = isfolder;
2033}
2034
2035////////////////////////////////////////////////////////////////////////////////
2036/// Set canvas name. In case `name` is an empty string, a default name is set.
2037
2038void TCanvas::SetName(const char *name)
2039{
2040 if (name && *name)
2041 fName = name;
2042 else
2044
2045 if (gPad && TestBit(kMustCleanup))
2046 gPad->Modified();
2047}
2048
2049
2050////////////////////////////////////////////////////////////////////////////////
2051/// Function to resize a canvas so that the plot inside is shown in real aspect
2052/// ratio
2053///
2054/// \param[in] axis 1 for resizing horizontally (x-axis) in order to get real
2055/// aspect ratio, 2 for the resizing vertically (y-axis)
2056/// \return false if error is encountered, true otherwise
2057///
2058/// ~~~ {.cpp}
2059/// hpxpy->Draw();
2060/// c1->SetRealAspectRatio();
2061/// ~~~
2062///
2063/// - For defining the concept of real aspect ratio, it is assumed that x and y
2064/// axes are in same units, e.g. both in MeV or both in ns.
2065/// - You can resize either the width of the canvas or the height, but not both
2066/// at the same time
2067/// - Call this function AFTER drawing AND zooming (SetUserRange) your TGraph or
2068/// Histogram, otherwise it cannot infer your actual axes lengths
2069/// - This function ensures that the TFrame has a real aspect ratio, this does not
2070/// mean that the full pad (i.e. the canvas or png output) including margins has
2071/// exactly the same ratio
2072/// - This function does not work if the canvas is divided in several subpads
2073
2074bool TCanvas::SetRealAspectRatio(const Int_t axis)
2075{
2076 Update();
2077
2078 //Get how many pixels are occupied by the canvas
2079 Int_t npx = GetWw();
2080 Int_t npy = GetWh();
2081
2082 //Get x-y coordinates at the edges of the canvas (extrapolating outside the axes, NOT at the edges of the histogram)
2083 Double_t x1 = GetX1();
2084 Double_t y1 = GetY1();
2085 Double_t x2 = GetX2();
2086 Double_t y2 = GetY2();
2087
2088 //Get the length of extrapolated x and y axes
2089 Double_t xlength2 = x2 - x1;
2090 Double_t ylength2 = y2 - y1;
2091 Double_t ratio2 = xlength2/ylength2;
2092
2093 //Now get the number of pixels including the canvas borders
2094 Int_t bnpx = GetWindowWidth();
2095 Int_t bnpy = GetWindowHeight();
2096
2097 if (axis==1) {
2098 SetCanvasSize(TMath::Nint(npy*ratio2), npy);
2099 SetWindowSize((bnpx-npx)+TMath::Nint(npy*ratio2), bnpy);
2100 } else if (axis==2) {
2101 SetCanvasSize(npx, TMath::Nint(npx/ratio2));
2102 SetWindowSize(bnpx, (bnpy-npy)+TMath::Nint(npx/ratio2));
2103 } else {
2104 Error("SetRealAspectRatio", "axis value %d is neither 1 (resize along x-axis) nor 2 (resize along y-axis).",axis);
2105 return false;
2106 }
2107
2108 //Check now that resizing has worked
2109
2110 Update();
2111
2112 //Get how many pixels are occupied by the canvas
2113 npx = GetWw();
2114 npy = GetWh();
2115
2116 //Get x-y coordinates at the edges of the canvas (extrapolating outside the axes,
2117 //NOT at the edges of the histogram)
2118 x1 = GetX1();
2119 y1 = GetY1();
2120 x2 = GetX2();
2121 y2 = GetY2();
2122
2123 //Get the length of extrapolated x and y axes
2124 xlength2 = x2 - x1;
2125 ylength2 = y2 - y1;
2126 ratio2 = xlength2/ylength2;
2127
2128 //Check accuracy +/-1 pixel due to rounding
2129 if (abs(TMath::Nint(npy*ratio2) - npx)<2) {
2130 return true;
2131 } else {
2132 Error("SetRealAspectRatio", "Resizing failed.");
2133 return false;
2134 }
2135}
2136
2137
2138////////////////////////////////////////////////////////////////////////////////
2139/// Set selected canvas.
2140
2142{
2143 fSelected = obj;
2144 if (obj) obj->SetBit(kMustCleanup);
2145}
2146
2147////////////////////////////////////////////////////////////////////////////////
2148/// Set canvas title.
2149
2150void TCanvas::SetTitle(const char *title)
2151{
2152 fTitle = title;
2154}
2155
2156////////////////////////////////////////////////////////////////////////////////
2157/// Set canvas window position
2158
2160{
2161 if (fCanvasImp)
2163}
2164
2165////////////////////////////////////////////////////////////////////////////////
2166/// Set canvas window size
2167
2169{
2170 if (fBatch)
2171 SetCanvasSize((ww + fCw) / 2, (wh + fCh) / 2);
2172 else if (fCanvasImp)
2173 fCanvasImp->SetWindowSize(ww, wh);
2174}
2175
2176////////////////////////////////////////////////////////////////////////////////
2177/// Set the canvas scale in centimeters.
2178///
2179/// This information is used by PostScript to set the page size.
2180///
2181/// \param[in] xsize size of the canvas in centimeters along X
2182/// \param[in] ysize size of the canvas in centimeters along Y
2183///
2184/// if xsize and ysize are not equal to 0, then the scale factors will
2185/// be computed to keep the ratio ysize/xsize independently of the canvas
2186/// size (parts of the physical canvas will be unused).
2187///
2188/// if xsize = 0 and ysize is not zero, then xsize will be computed
2189/// to fit to the current canvas scale. If the canvas is resized,
2190/// a new value for xsize will be recomputed. In this case the aspect
2191/// ratio is not preserved.
2192///
2193/// if both xsize = 0 and ysize = 0, then the scaling is automatic.
2194/// the largest dimension will be allocated a size of 20 centimeters.
2195
2196void TCanvas::Size(Float_t xsize, Float_t ysize)
2197{
2198 fXsizeUser = xsize;
2199 fYsizeUser = ysize;
2200
2201 Resize();
2202}
2203
2204////////////////////////////////////////////////////////////////////////////////
2205/// Show canvas
2206
2207void TCanvas::Show()
2208{
2209 if (fCanvasImp)
2210 fCanvasImp->Show();
2211}
2212
2213////////////////////////////////////////////////////////////////////////////////
2214/// Stream a class object.
2215
2217{
2218 UInt_t R__s, R__c;
2219 if (b.IsReading()) {
2220 Version_t v = b.ReadVersion(&R__s, &R__c);
2221 gPad = this;
2222 fCanvas = this;
2223 if (v>7) b.ClassBegin(TCanvas::IsA());
2224 if (v>7) b.ClassMember("TPad");
2226 gPad = this;
2227 //restore the colors
2228 TObjArray *colors = (TObjArray*)fPrimitives->FindObject("ListOfColors");
2229 if (colors) {
2230 TIter next(colors);
2231 TColor *colold;
2232 while ((colold = (TColor*)next())) {
2233 if (colold) {
2234 Int_t cn = 0;
2235 if (colold) cn = colold->GetNumber();
2236 TColor *colcur = gROOT->GetColor(cn);
2237 if (colcur) {
2238 colcur->SetRGB(colold->GetRed(),colold->GetGreen(),colold->GetBlue());
2239 } else {
2240 colcur = new TColor(cn,colold->GetRed(),
2241 colold->GetGreen(),
2242 colold->GetBlue(),
2243 colold->GetName());
2244 if (!colcur) return;
2245 }
2246 }
2247 }
2248 //restore the palette if needed
2249 TObjArray *currentpalette = (TObjArray*)fPrimitives->FindObject("CurrentColorPalette");
2250 if (currentpalette) {
2251 TIter nextpal(currentpalette);
2252 Int_t n = currentpalette->GetEntries();
2253 TArrayI palcolors(n);
2254 TColor *col = nullptr;
2255 Int_t i = 0;
2256 while ((col = (TColor*)nextpal())) palcolors[i++] = col->GetNumber();
2257 gStyle->SetPalette(n,palcolors.GetArray());
2258 fPrimitives->Remove(currentpalette);
2259 delete currentpalette;
2260 }
2262 colors->Delete();
2263 delete colors;
2264 }
2265
2266 if (v>7) b.ClassMember("fDISPLAY","TString");
2268 if (v>7) b.ClassMember("fDoubleBuffer", "Int_t");
2269 b >> fDoubleBuffer;
2270 if (v>7) b.ClassMember("fRetained", "Bool_t");
2271 b >> fRetained;
2272 if (v>7) b.ClassMember("fXsizeUser", "Size_t");
2273 b >> fXsizeUser;
2274 if (v>7) b.ClassMember("fYsizeUser", "Size_t");
2275 b >> fYsizeUser;
2276 if (v>7) b.ClassMember("fXsizeReal", "Size_t");
2277 b >> fXsizeReal;
2278 if (v>7) b.ClassMember("fYsizeReal", "Size_t");
2279 b >> fYsizeReal;
2280 fCanvasID = -1;
2281 if (v>7) b.ClassMember("fWindowTopX", "Int_t");
2282 b >> fWindowTopX;
2283 if (v>7) b.ClassMember("fWindowTopY", "Int_t");
2284 b >> fWindowTopY;
2285 if (v > 2) {
2286 if (v>7) b.ClassMember("fWindowWidth", "UInt_t");
2287 b >> fWindowWidth;
2288 if (v>7) b.ClassMember("fWindowHeight", "UInt_t");
2289 b >> fWindowHeight;
2290 }
2291 if (v>7) b.ClassMember("fCw", "UInt_t");
2292 b >> fCw;
2293 if (v>7) b.ClassMember("fCh", "UInt_t");
2294 b >> fCh;
2295 if (v <= 2) {
2296 fWindowWidth = fCw;
2298 }
2299 if (v>7) b.ClassMember("fCatt", "TAttCanvas");
2300 fCatt.Streamer(b);
2301 Bool_t dummy;
2302 if (v>7) b.ClassMember("kMoveOpaque", "Bool_t");
2303 b >> dummy; if (dummy) MoveOpaque(1);
2304 if (v>7) b.ClassMember("kResizeOpaque", "Bool_t");
2305 b >> dummy; if (dummy) ResizeOpaque(1);
2306 if (v>7) b.ClassMember("fHighLightColor", "Color_t");
2307 b >> fHighLightColor;
2308 if (v>7) b.ClassMember("fBatch", "Bool_t");
2309 b >> dummy; //was fBatch
2310 if (v < 2) return;
2311 if (v>7) b.ClassMember("kShowEventStatus", "Bool_t");
2312 b >> dummy; if (dummy) SetBit(kShowEventStatus);
2313
2314 if (v > 3) {
2315 if (v>7) b.ClassMember("kAutoExec", "Bool_t");
2316 b >> dummy; if (dummy) SetBit(kAutoExec);
2317 }
2318 if (v>7) b.ClassMember("kMenuBar", "Bool_t");
2319 b >> dummy; if (dummy) SetBit(kMenuBar);
2320 fBatch = gROOT->IsBatch();
2321 if (v>7) b.ClassEnd(TCanvas::IsA());
2322 b.CheckByteCount(R__s, R__c, TCanvas::IsA());
2323 } else {
2324 //save list of colors
2325 //we must protect the case when two or more canvases are saved
2326 //in the same buffer. If the list of colors has already been saved
2327 //in the buffer, do not add the list of colors to the list of primitives.
2328 TObjArray *colors = nullptr;
2329 TObjArray *CurrentColorPalette = nullptr;
2330 if (TColor::DefinedColors()) {
2331 if (!b.CheckObject(gROOT->GetListOfColors(),TObjArray::Class())) {
2332 colors = (TObjArray*)gROOT->GetListOfColors();
2334 }
2335 //save the current palette
2337 Int_t palsize = pal.GetSize();
2338 CurrentColorPalette = new TObjArray();
2339 CurrentColorPalette->SetName("CurrentColorPalette");
2340 for (Int_t i=0; i<palsize; i++) CurrentColorPalette->Add(gROOT->GetColor(pal[i]));
2341 fPrimitives->Add(CurrentColorPalette);
2342 }
2343
2344 R__c = b.WriteVersion(TCanvas::IsA(), kTRUE);
2345 b.ClassBegin(TCanvas::IsA());
2346 b.ClassMember("TPad");
2349 if (CurrentColorPalette) { fPrimitives->Remove(CurrentColorPalette); delete CurrentColorPalette; }
2350 b.ClassMember("fDISPLAY","TString");
2352 b.ClassMember("fDoubleBuffer", "Int_t");
2353 b << fDoubleBuffer;
2354 b.ClassMember("fRetained", "Bool_t");
2355 b << fRetained;
2356 b.ClassMember("fXsizeUser", "Size_t");
2357 b << fXsizeUser;
2358 b.ClassMember("fYsizeUser", "Size_t");
2359 b << fYsizeUser;
2360 b.ClassMember("fXsizeReal", "Size_t");
2361 b << fXsizeReal;
2362 b.ClassMember("fYsizeReal", "Size_t");
2363 b << fYsizeReal;
2365 Int_t topx = fWindowTopX, topy = fWindowTopY;
2366 UInt_t editorWidth = 0;
2367 if(fCanvasImp) editorWidth = fCanvasImp->GetWindowGeometry(topx,topy,w,h);
2368 b.ClassMember("fWindowTopX", "Int_t");
2369 b << topx;
2370 b.ClassMember("fWindowTopY", "Int_t");
2371 b << topy;
2372 b.ClassMember("fWindowWidth", "UInt_t");
2373 b << (UInt_t)(w-editorWidth);
2374 b.ClassMember("fWindowHeight", "UInt_t");
2375 b << h;
2376 b.ClassMember("fCw", "UInt_t");
2377 b << fCw;
2378 b.ClassMember("fCh", "UInt_t");
2379 b << fCh;
2380 b.ClassMember("fCatt", "TAttCanvas");
2381 fCatt.Streamer(b);
2382 b.ClassMember("kMoveOpaque", "Bool_t");
2383 b << TestBit(kMoveOpaque); //please remove in ROOT version 6
2384 b.ClassMember("kResizeOpaque", "Bool_t");
2385 b << TestBit(kResizeOpaque); //please remove in ROOT version 6
2386 b.ClassMember("fHighLightColor", "Color_t");
2387 b << fHighLightColor;
2388 b.ClassMember("fBatch", "Bool_t");
2389 b << fBatch; //please remove in ROOT version 6
2390 b.ClassMember("kShowEventStatus", "Bool_t");
2391 b << TestBit(kShowEventStatus); //please remove in ROOT version 6
2392 b.ClassMember("kAutoExec", "Bool_t");
2393 b << TestBit(kAutoExec); //please remove in ROOT version 6
2394 b.ClassMember("kMenuBar", "Bool_t");
2395 b << TestBit(kMenuBar); //please remove in ROOT version 6
2396 b.ClassEnd(TCanvas::IsA());
2397 b.SetByteCount(R__c, kTRUE);
2398 }
2399}
2400
2401////////////////////////////////////////////////////////////////////////////////
2402/// Toggle pad auto execution of list of TExecs.
2403
2405{
2406 Bool_t autoExec = TestBit(kAutoExec);
2407 SetBit(kAutoExec,!autoExec);
2408}
2409
2410////////////////////////////////////////////////////////////////////////////////
2411/// Toggle event statusbar.
2412
2414{
2415 Bool_t showEventStatus = !TestBit(kShowEventStatus);
2416 SetBit(kShowEventStatus,showEventStatus);
2417
2418 if (fCanvasImp) fCanvasImp->ShowStatusBar(showEventStatus);
2419}
2420
2421////////////////////////////////////////////////////////////////////////////////
2422/// Toggle toolbar.
2423
2425{
2426 Bool_t showToolBar = !TestBit(kShowToolBar);
2427 SetBit(kShowToolBar,showToolBar);
2428
2429 if (fCanvasImp) fCanvasImp->ShowToolBar(showToolBar);
2430}
2431
2432////////////////////////////////////////////////////////////////////////////////
2433/// Toggle editor.
2434
2436{
2437 Bool_t showEditor = !TestBit(kShowEditor);
2438 SetBit(kShowEditor,showEditor);
2439
2440 if (fCanvasImp) fCanvasImp->ShowEditor(showEditor);
2441}
2442
2443////////////////////////////////////////////////////////////////////////////////
2444/// Toggle tooltip display.
2445
2447{
2448 Bool_t showToolTips = !TestBit(kShowToolTips);
2449 SetBit(kShowToolTips, showToolTips);
2450
2451 if (fCanvasImp) fCanvasImp->ShowToolTips(showToolTips);
2452}
2453
2454
2455////////////////////////////////////////////////////////////////////////////////
2456/// Static function returning "true" if transparency is supported.
2457
2459{
2460 return gPad && (gVirtualX->InheritsFrom("TGQuartz") ||
2461 gPad->GetGLDevice() != -1);
2462}
2463
2464extern "C" void ROOT_TCanvas_Update(void* TheCanvas) {
2465 static_cast<TCanvas*>(TheCanvas)->Update();
2466}
2467
2468////////////////////////////////////////////////////////////////////////////////
2469/// Update canvas pad buffers.
2470
2471void TCanvas::Update()
2472{
2473 if (fUpdating) return;
2474
2475 if (fPixmapID == -1) return;
2476
2477 static const union CastFromFuncToVoidPtr_t {
2478 CastFromFuncToVoidPtr_t(): fFuncPtr(ROOT_TCanvas_Update) {}
2479 void (*fFuncPtr)(void*);
2480 void* fVoidPtr;
2481 } castFromFuncToVoidPtr;
2482
2483 if (gThreadXAR) {
2484 void *arr[3];
2485 arr[1] = this;
2486 arr[2] = castFromFuncToVoidPtr.fVoidPtr;
2487 if ((*gThreadXAR)("CUPD", 3, arr, nullptr)) return;
2488 }
2489
2490 if (!fCanvasImp) return;
2491
2492 if (!gVirtualX->IsCmdThread()) {
2493 // Why do we have this (which uses the interpreter to funnel the Update()
2494 // through the main thread) when the gThreadXAR mechanism does seemingly
2495 // the same?
2496 gInterpreter->Execute(this, IsA(), "Update", "");
2497 return;
2498 }
2499
2501
2502 fUpdating = kTRUE;
2503
2504 if (!fCanvasImp->PerformUpdate()) {
2505
2506 if (!IsBatch()) FeedbackMode(kFALSE); // Goto double buffer mode
2507
2508 if (!UseGL() || fGLDevice == -1) PaintModified(); // Repaint all modified pad's
2509
2510 Flush(); // Copy all pad pixmaps to the screen
2511
2513 }
2514
2515 fUpdating = kFALSE;
2516}
2517
2518////////////////////////////////////////////////////////////////////////////////
2519/// Used by friend class TCanvasImp.
2520
2522{
2523 fCanvasID = 0;
2524 fContextMenu = nullptr;
2525}
2526
2527////////////////////////////////////////////////////////////////////////////////
2528/// Check whether this canvas is to be drawn in grayscale mode.
2529
2531{
2532 return TestBit(kIsGrayscale);
2533}
2534
2535////////////////////////////////////////////////////////////////////////////////
2536/// Set whether this canvas should be painted in grayscale, and re-paint
2537/// it if necessary.
2538
2539void TCanvas::SetGrayscale(Bool_t set /*= kTRUE*/)
2540{
2541 if (IsGrayscale() == set) return;
2542 SetBit(kIsGrayscale, set);
2543 Paint(); // update canvas and all sub-pads, unconditionally!
2544}
2545
2546////////////////////////////////////////////////////////////////////////////////
2547/// Probably, TPadPainter must be placed in a separate ROOT module -
2548/// "padpainter" (the same as "histpainter"). But now, it's directly in a
2549/// gpad dir, so, in case of default painter, no *.so should be loaded,
2550/// no need in plugin managers.
2551/// May change in future.
2552
2554{
2555 //Even for batch mode painter is still required, just to delegate
2556 //some calls to batch "virtual X".
2557 if (!UseGL() || fBatch) {
2558 fPainter = nullptr;
2560 if (!fPainter) fPainter = new TPadPainter; // Do not need plugin manager for this!
2561 } else {
2563 if (!fPainter) {
2564 Error("CreatePainter", "GL Painter creation failed! Will use default!");
2565 fPainter = new TPadPainter;
2566 fUseGL = kFALSE;
2567 }
2568 }
2569}
2570
2571////////////////////////////////////////////////////////////////////////////////
2572/// Access and (probably) creation of pad painter.
2573
2575{
2576 if (!fPainter) CreatePainter();
2577 return fPainter;
2578}
2579
2580
2581////////////////////////////////////////////////////////////////////////////////
2582///assert on IsBatch() == false?
2583
2585{
2586 if (fGLDevice != -1) {
2587 //fPainter has a font manager.
2588 //Font manager will delete textures.
2589 //If context is wrong (we can have several canvases) -
2590 //wrong texture will be deleted, damaging some of our fonts.
2591 gGLManager->MakeCurrent(fGLDevice);
2592 }
2593
2595
2596 if (fGLDevice != -1) {
2597 gGLManager->DeleteGLContext(fGLDevice);//?
2598 fGLDevice = -1;
2599 }
2600}
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
ECursor
Definition GuiTypes.h:372
@ kCross
Definition GuiTypes.h:374
#define SafeDelete(p)
Definition RConfig.hxx:540
#define b(i)
Definition RSha256.hxx:100
#define c(i)
Definition RSha256.hxx:101
#define h(i)
Definition RSha256.hxx:106
bool Bool_t
Definition RtypesCore.h:63
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
float Float_t
Definition RtypesCore.h:57
constexpr Bool_t kFALSE
Definition RtypesCore.h:101
double Double_t
Definition RtypesCore.h:59
constexpr 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:2461
class TCanvasInit gCanvasInit
auto GetNewCanvasName()
Definition TCanvas.cxx:66
const Size_t kDefaultCanvasSize
Definition TCanvas.cxx:61
R__EXTERN TEnv * gEnv
Definition TEnv.h:170
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t option
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t cursor
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void SetDoubleBuffer
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char filename
Option_t Option_t SetFillStyle
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char cname
Option_t Option_t TPoint TPoint const char x2
Option_t Option_t TPoint TPoint const char x1
Option_t Option_t TPoint TPoint const char mode
Option_t Option_t TPoint TPoint const char y2
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void SetCursor
Option_t Option_t SetFillColor
Option_t Option_t TPoint TPoint const char y1
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:283
Int_t gDebug
Definition TROOT.cxx:585
R__EXTERN TVirtualMutex * gROOTMutex
Definition TROOT.h:63
#define gROOT
Definition TROOT.h:405
R__EXTERN TStyle * gStyle
Definition TStyle.h:414
#define gGLManager
Definition TVirtualGL.h:160
#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(Bool_t only_web=kFALSE)
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
const Int_t * GetArray() const
Definition TArrayI.h:43
Int_t GetSize() const
Definition TArray.h:47
virtual void Streamer(TBuffer &)
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
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
static TClass * Class()
virtual Bool_t GetTimeDisplay() const
Definition TAxis.h:127
virtual const char * GetTimeFormat() const
Definition TAxis.h:128
Create a Box.
Definition TBox.h:22
static TClass * Class()
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 SetStatusText(const char *text=nullptr, Int_t partidx=0)
Definition TCanvasImp.h:69
virtual void SetWindowPosition(Int_t x, Int_t y)
Definition TCanvasImp.h:70
virtual void ShowMenuBar(Bool_t show=kTRUE)
Definition TCanvasImp.h:75
virtual void Show()
Definition TCanvasImp.h:74
virtual void ShowToolTips(Bool_t show=kTRUE)
Definition TCanvasImp.h:82
virtual void Iconify()
Definition TCanvasImp.h:67
virtual Int_t InitWindow()
Definition TCanvasImp.h:68
virtual void Close()
Definition TCanvasImp.h:59
virtual void SetWindowTitle(const char *newTitle)
Definition TCanvasImp.h:72
virtual Bool_t IsWeb() const
Definition TCanvasImp.h:48
virtual UInt_t GetWindowGeometry(Int_t &x, Int_t &y, UInt_t &w, UInt_t &h)
Definition TCanvasImp.h:61
virtual TVirtualPadPainter * CreatePadPainter()
Definition TCanvasImp.h:50
virtual void ShowEditor(Bool_t show=kTRUE)
Definition TCanvasImp.h:80
virtual void RaiseWindow()
Definition TCanvasImp.h:77
virtual Bool_t PerformUpdate()
Definition TCanvasImp.h:49
virtual void SetWindowSize(UInt_t width, UInt_t height)
Definition TCanvasImp.h:71
virtual void SetCanvasSize(UInt_t w, UInt_t h)
Definition TCanvasImp.h:73
virtual void ForceUpdate()
Definition TCanvasImp.h:60
virtual void ShowStatusBar(Bool_t show=kTRUE)
Definition TCanvasImp.h:76
virtual void ShowToolBar(Bool_t show=kTRUE)
Definition TCanvasImp.h:81
The Canvas class.
Definition TCanvas.h:23
void Init()
Initialize the TCanvas members. Called by all constructors.
Definition TCanvas.cxx:533
UInt_t fCw
Width of the canvas along X (pixels)
Definition TCanvas.h:43
void EmbedInto(Int_t winid, Int_t ww, Int_t wh)
Embedded a canvas into a TRootEmbeddedCanvas.
Definition TCanvas.cxx:1040
void SetWindowSize(UInt_t ww, UInt_t wh)
Set canvas window size.
Definition TCanvas.cxx:2165
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:2027
void Browse(TBrowser *b) override
Browse.
Definition TCanvas.cxx:680
UInt_t GetWindowHeight() const
Definition TCanvas.h:161
virtual void EditorBar()
Get editor bar.
Definition TCanvas.cxx:1031
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:1066
Size_t fYsizeReal
Current size of canvas along Y in CM.
Definition TCanvas.h:36
void Constructor()
Canvas default constructor.
Definition TCanvas.cxx:186
virtual void ToggleAutoExec()
Toggle pad auto execution of list of TExecs.
Definition TCanvas.cxx:2401
TCanvas(const TCanvas &canvas)=delete
Int_t fWindowTopX
Top X position of window (in pixels)
Definition TCanvas.h:39
void Draw(Option_t *option="") override
Draw a canvas.
Definition TCanvas.cxx:850
void SetDoubleBuffer(Int_t mode=1) override
Set Double Buffer On/Off.
Definition TCanvas.cxx:1988
virtual void ToggleToolTips()
Toggle tooltip display.
Definition TCanvas.cxx:2443
void Clear(Option_t *option="") override
Remove all primitives from the canvas.
Definition TCanvas.cxx:733
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:2410
void Destructor()
Actual canvas destructor.
Definition TCanvas.cxx:690
void DeleteCanvasPainter()
assert on IsBatch() == false?
Definition TCanvas.cxx:2581
TPad * fPadSave
! Pointer to saved pad in HandleInput
Definition TCanvas.h:56
static Bool_t SupportAlpha()
Static function returning "true" if transparency is supported.
Definition TCanvas.cxx:2455
Bool_t fBatch
! True when in batchmode
Definition TCanvas.h:59
Bool_t fUseGL
! True when rendering is with GL
Definition TCanvas.h:62
Int_t fEventX
! Last X mouse position in canvas
Definition TCanvas.h:46
Bool_t IsBatch() const override
Definition TCanvas.h:170
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:898
Size_t fXsizeReal
Current size of canvas along X in CM.
Definition TCanvas.h:35
Bool_t HasMenuBar() const
Definition TCanvas.h:167
TVirtualPadPainter * GetCanvasPainter()
Access and (probably) creation of pad painter.
Definition TCanvas.cxx:2571
virtual void HighlightConnect(const char *slot)
This is "simplification" for function TCanvas::Connect with Highlighted signal for specific slot.
Definition TCanvas.cxx:1619
TPad * Pick(Int_t px, Int_t py, TObjLink *&pickobj) override
Search for an object at pixel position px,py.
Definition TCanvas.h:181
void Close(Option_t *option="") override
Close canvas.
Definition TCanvas.cxx:784
void SetFixedAspectRatio(Bool_t fixed=kTRUE) override
Fix canvas aspect ratio to current value if fixed is true.
Definition TCanvas.cxx:2005
virtual void Resize(Option_t *option="")
Recompute canvas parameters following a X11 Resize.
Definition TCanvas.cxx:1656
Color_t GetHighLightColor() const override
Definition TCanvas.h:137
Bool_t GetShowToolBar() const
Definition TCanvas.h:148
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:970
Bool_t IsFolder() const override
Is folder ?
Definition TCanvas.cxx:1476
UInt_t fWindowWidth
Width of window (including borders, etc.)
Definition TCanvas.h:41
TVirtualPadPainter * fPainter
! Canvas (pad) painter.
Definition TCanvas.h:65
void CopyPixmaps() override
Copy the canvas pixmap of the pad to the canvas.
Definition TCanvas.cxx:832
Bool_t IsGrayscale()
Check whether this canvas is to be drawn in grayscale mode.
Definition TCanvas.cxx:2527
TPad * fClickSelectedPad
! Pad containing currently click-selected object
Definition TCanvas.h:55
Bool_t fUpdating
! True when Updating the canvas
Definition TCanvas.h:60
void SaveSource(const char *filename="", Option_t *option="")
Save primitives in this canvas as a C++ macro file.
Definition TCanvas.cxx:1810
Color_t fHighLightColor
Highlight color of active pad.
Definition TCanvas.h:37
virtual void Size(Float_t xsizeuser=0, Float_t ysizeuser=0)
Set the canvas scale in centimeters.
Definition TCanvas.cxx:2193
virtual void ProcessedEvent(Int_t event, Int_t x, Int_t y, TObject *selected)
Emit ProcessedEvent() signal.
Definition TCanvas.cxx:1641
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:33
Int_t fEventY
! Last Y mouse position in canvas
Definition TCanvas.h:47
TVirtualPad * cd(Int_t subpadnumber=0) override
Set current canvas & pad.
Definition TCanvas.cxx:714
UInt_t fWindowHeight
Height of window (including menubar, borders, etc.)
Definition TCanvas.h:42
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:50
void SetCanvasSize(UInt_t ww, UInt_t wh) override
Set Width and Height of canvas to ww and wh respectively.
Definition TCanvas.cxx:1965
void Show()
Show canvas.
Definition TCanvas.cxx:2204
TPad * fSelectedPad
! Pad containing currently selected object
Definition TCanvas.h:54
virtual void Selected(TVirtualPad *pad, TObject *obj, Int_t event)
Emit Selected() signal.
Definition TCanvas.cxx:1627
Int_t fSelectedX
! X of selected object
Definition TCanvas.h:51
virtual void ToggleEditor()
Toggle editor.
Definition TCanvas.cxx:2432
TVirtualPad * GetSelectedPad() const override
Definition TCanvas.h:145
virtual void Picked(TPad *selpad, TObject *selected, Int_t event)
Emit Picked() signal.
Definition TCanvas.cxx:1581
TObject * fSelected
! Currently selected object
Definition TCanvas.h:49
void SavePrimitive(std::ostream &out, Option_t *option="") override
Save primitives in this canvas in C++ macro file with GUI.
Definition TCanvas.cxx:1766
void SetCursor(ECursor cursor) override
Set cursor.
Definition TCanvas.cxx:1979
void Streamer(TBuffer &) override
Stream a class object.
Definition TCanvas.cxx:2213
Bool_t GetShowToolTips() const
Definition TCanvas.h:150
Int_t fCanvasID
! Canvas identifier
Definition TCanvas.h:48
void SetGrayscale(Bool_t set=kTRUE)
Set whether this canvas should be painted in grayscale, and re-paint it if necessary.
Definition TCanvas.cxx:2536
void SetTitle(const char *title="") override
Set canvas title.
Definition TCanvas.cxx:2147
UInt_t fCh
Height of the canvas along Y (pixels)
Definition TCanvas.h:44
TContextMenu * fContextMenu
! Context menu pointer
Definition TCanvas.h:58
TAttCanvas fCatt
Canvas attributes.
Definition TCanvas.h:31
void SetName(const char *name="") override
Set canvas name. In case name is an empty string, a default name is set.
Definition TCanvas.cxx:2035
UInt_t GetWindowWidth() const
Definition TCanvas.h:160
Bool_t fRetained
Retain structure flag.
Definition TCanvas.h:61
void DisconnectWidget()
Used by friend class TCanvasImp.
Definition TCanvas.cxx:2518
void FeedbackMode(Bool_t set)
Turn rubberband feedback mode on or off.
Definition TCanvas.cxx:1120
void ls(Option_t *option="") const override
List all pads.
Definition TCanvas.cxx:1492
void RaiseWindow()
Raise canvas window.
Definition TCanvas.cxx:1732
void Build()
Build a canvas. Called by all constructors.
Definition TCanvas.cxx:582
Int_t fWindowTopY
Top Y position of window (in pixels)
Definition TCanvas.h:40
void Paint(Option_t *option="") override
Paint canvas.
Definition TCanvas.cxx:1531
TClass * IsA() const override
Definition TCanvas.h:232
void Update() override
Update canvas pad buffers.
Definition TCanvas.cxx:2468
void ExecuteEvent(Int_t event, Int_t px, Int_t py) override
Execute action corresponding to one event.
Definition TCanvas.cxx:1102
void RunAutoExec()
Execute the list of TExecs in the current pad.
Definition TCanvas.cxx:1755
virtual void Cleared(TVirtualPad *pad)
Emit pad Cleared signal.
Definition TCanvas.cxx:766
UInt_t GetWw() const override
Definition TCanvas.h:162
TCanvasImp * fCanvasImp
! Window system specific canvas implementation
Definition TCanvas.h:57
UInt_t GetWh() const override
Definition TCanvas.h:163
virtual void Highlighted(TVirtualPad *pad, TObject *obj, Int_t x, Int_t y)
Emit Highlighted() signal.
Definition TCanvas.cxx:1600
void Flush()
Flush canvas buffers.
Definition TCanvas.cxx:1134
Size_t fYsizeUser
User specified size of canvas along Y in CM.
Definition TCanvas.h:34
Int_t fDoubleBuffer
Double buffer flag (0=off, 1=on)
Definition TCanvas.h:38
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:2550
void SetSelected(TObject *obj) override
Set selected canvas.
Definition TCanvas.cxx:2138
void MoveOpaque(Int_t set=1)
Set option to move objects/pads in a canvas.
Definition TCanvas.cxx:1523
virtual ~TCanvas()
Canvas destructor.
Definition TCanvas.cxx:672
static Bool_t fgIsFolder
Indicates if canvas can be browsed as a folder.
Definition TCanvas.h:67
void Closed() override
Emit Closed signal.
Definition TCanvas.cxx:774
void SetWindowPosition(Int_t x, Int_t y)
Set canvas window position.
Definition TCanvas.cxx:2156
TString fDISPLAY
Name of destination screen.
Definition TCanvas.h:32
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:2071
Int_t fEvent
! Type of current or last handled event
Definition TCanvas.h:45
Bool_t GetShowEventStatus() const
Definition TCanvas.h:147
TString fSelectedOpt
! Drawing option of selected object
Definition TCanvas.h:53
Int_t fSelectedY
! Y of selected object
Definition TCanvas.h:52
Bool_t fDrawn
! Set to True when the Draw method is called
Definition TCanvas.h:63
void SetBatch(Bool_t batch=kTRUE) override
Toggle batch mode.
Definition TCanvas.cxx:1949
Bool_t UseGL() const
Definition TCanvas.h:224
@ kResizeOpaque
Definition TCanvas.h:94
@ kShowToolTips
Definition TCanvas.h:96
@ kShowToolBar
Definition TCanvas.h:91
@ kMoveOpaque
Definition TCanvas.h:93
@ kIsGrayscale
Definition TCanvas.h:95
@ kShowEventStatus
Definition TCanvas.h:88
@ kAutoExec
Definition TCanvas.h:89
@ kMenuBar
Definition TCanvas.h:90
@ kShowEditor
Definition TCanvas.h:92
void ResizeOpaque(Int_t set=1)
Set option to resize objects/pads in a canvas.
Definition TCanvas.cxx:1747
virtual void ToggleToolBar()
Toggle toolbar.
Definition TCanvas.cxx:2421
virtual TObject * DrawClonePad()
Draw a clone of this canvas into the current pad In an interactive session, select the destination/cu...
Definition TCanvas.cxx:915
TClass instances represent classes, structs and namespaces in the ROOT type system.
Definition TClass.h:81
@ kRealNew
Definition TClass.h:108
static ENewType IsCallingNew()
Static method returning the defConstructor flag passed to TClass::New().
Definition TClass.cxx:5902
@ kClassSaved
Definition TClass.h:95
void SetName(const char *name)
void Browse(TBrowser *b) override
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:1759
static const TArrayI & GetPalette()
Static function returning the current active palette.
Definition TColor.cxx:1462
Float_t GetRed() const
Definition TColor.h:58
static Bool_t 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:2188
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:1823
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:1480
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
static TClass * Class()
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
TObject * FindObject(const char *name) const override
Find an object in this list using its name.
Definition TList.cxx:578
void Add(TObject *obj) override
Definition TList.h:81
TObject * Remove(TObject *obj) override
Remove object from the list.
Definition TList.cxx:822
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
An array of TObjects.
Definition TObjArray.h:31
Int_t GetEntries() const override
Return the number of objects in array (i.e.
static TClass * Class()
void Add(TObject *obj) override
Definition TObjArray.h:68
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:439
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:223
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:207
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:956
virtual void ExecuteEvent(Int_t event, Int_t px, Int_t py)
Execute action corresponding to an event at (px,py).
Definition TObject.cxx:393
virtual void Execute(const char *method, const char *params, Int_t *error=nullptr)
Execute method on this object with the given parameter string, e.g.
Definition TObject.cxx:359
virtual char * GetObjectInfo(Int_t px, Int_t py) const
Returns string containing info about the object at position (px,py).
Definition TObject.cxx:468
virtual void Delete(Option_t *option="")
Delete this object.
Definition TObject.cxx:248
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:774
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition TObject.cxx:525
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:970
virtual const char * GetTitle() const
Returns title of object.
Definition TObject.cxx:483
virtual void Pop()
Pop on object drawn in a pad to the top of the display list.
Definition TObject.cxx:610
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:944
Implement TVirtualPadPainter which abstracts painting operations.
Definition TPadPainter.h:26
The most important graphics class in the ROOT system.
Definition TPad.h:28
Short_t GetBorderMode() const override
Definition TPad.h:197
void SetBorderSize(Short_t bordersize) override
Definition TPad.h:321
Int_t GetTicky() const override
Definition TPad.h:237
void PaintBorder(Color_t color, Bool_t tops)
Paint the pad border.
Definition TPad.cxx:3508
void ResizePad(Option_t *option="") override
Compute pad conversion coefficients.
Definition TPad.cxx:5471
void SetGrid(Int_t valuex=1, Int_t valuey=1) override
Definition TPad.h:330
void SavePrimitive(std::ostream &out, Option_t *option="") override
Save primitives in this pad on the C++ source file out.
Definition TPad.cxx:5694
Bool_t GetGridx() const override
Definition TPad.h:233
Double_t fX2
X of upper X coordinate.
Definition TPad.h:38
void SetLogz(Int_t value=1) override
Set Lin/Log scale for Z.
Definition TPad.cxx:5945
Double_t GetY2() const override
Definition TPad.h:241
void Close(Option_t *option="") override
Delete all primitives in pad and pad itself.
Definition TPad.cxx:974
void PaintModified() override
Traverse pad hierarchy and (re)paint only modified pads.
Definition TPad.cxx:3676
const char * GetTitle() const override
Returns title of object.
Definition TPad.h:259
Double_t fX1
X of lower X coordinate.
Definition TPad.h:36
TList * GetListOfPrimitives() const override
Definition TPad.h:243
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:6005
void UseCurrentStyle() override
Force a copy of current style for all objects in pad.
Definition TPad.cxx:6758
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:5193
Double_t fY1
Y of lower Y coordinate.
Definition TPad.h:37
Int_t fGLDevice
! OpenGL off-screen pixmap identifier
Definition TPad.h:85
void Clear(Option_t *option="") override
Delete all pad primitives.
Definition TPad.cxx:626
Int_t GetTickx() const override
Definition TPad.h:236
Double_t fAspectRatio
ratio of w/h in case of fixed ratio
Definition TPad.h:82
void Modified(Bool_t flag=1) override
Definition TPad.h:417
void SetLogy(Int_t value=1) override
Set Lin/Log scale for Y.
Definition TPad.cxx:5934
TCanvas * fCanvas
! Pointer to mother canvas
Definition TPad.h:106
Bool_t fFixedAspectRatio
True if fixed aspect ratio.
Definition TPad.h:104
void ls(Option_t *option="") const override
List all primitives in pad.
Definition TPad.cxx:2966
void Streamer(TBuffer &) override
Stream a class object.
Definition TPad.cxx:6547
TString fTitle
Pad title.
Definition TPad.h:110
void CopyPixmaps() override
Copy the sub-pixmaps of the pad to the canvas.
Definition TPad.cxx:1047
void CopyPixmap() override
Copy the pixmap of the pad to the canvas.
Definition TPad.cxx:1033
Double_t GetY1() const override
Definition TPad.h:240
Bool_t GetGridy() const override
Definition TPad.h:234
void ExecuteEvent(Int_t event, Int_t px, Int_t py) override
Execute action corresponding to one event.
Definition TPad.cxx:1685
Int_t GetLogz() const override
Definition TPad.h:256
Short_t GetBorderSize() const override
Definition TPad.h:198
TList * fPrimitives
->List of primitives (subpads)
Definition TPad.h:107
TCanvas * GetCanvas() const override
Definition TPad.h:260
Short_t fBorderSize
pad bordersize in pixels
Definition TPad.h:97
void Paint(Option_t *option="") override
Paint all primitives in pad.
Definition TPad.cxx:3448
TString fName
Pad name.
Definition TPad.h:109
Int_t fPixmapID
! Off-screen pixmap identifier
Definition TPad.h:84
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:2595
TVirtualPad * cd(Int_t subpadnumber=0) override
Set Current pad.
Definition TPad.cxx:597
Int_t GetLogy() const override
Definition TPad.h:255
void SetBorderMode(Short_t bordermode) override
Definition TPad.h:320
void SetTicks(Int_t valuex=1, Int_t valuey=1) override
Definition TPad.h:350
Double_t fY2
Y of upper Y coordinate.
Definition TPad.h:39
Short_t fBorderMode
Bordermode (-1=down, 0 = no border, 1=up)
Definition TPad.h:98
void SetLogx(Int_t value=1) override
Set Lin/Log scale for X.
Definition TPad.cxx:5920
Int_t GetLogx() const override
Definition TPad.h:254
Double_t GetX2() const override
Definition TPad.h:239
Double_t GetX1() const override
Definition TPad.h:238
TPad * fMother
! pointer to mother of the list
Definition TPad.h:105
const char * GetName() const override
Returns name of object.
Definition TPad.h:258
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:2826
static void IndentLevel()
Functions used by ls() to indent an object hierarchy.
Definition TROOT.cxx:2834
static Int_t DecreaseDirLevel()
Decrease the indentation level for ls().
Definition TROOT.cxx:2713
Basic string class.
Definition TString.h:139
Ssiz_t Length() const
Definition TString.h:421
void ToLower()
Change string to lower-case.
Definition TString.cxx:1170
TSubString Strip(EStripType s=kTrailing, char c=' ') const
Return a substring of self stripped at beginning and/or end.
Definition TString.cxx:1151
const char * Data() const
Definition TString.h:380
@ kBoth
Definition TString.h:278
Ssiz_t Last(char c) const
Find last occurrence of a character c.
Definition TString.cxx:924
Bool_t IsNull() const
Definition TString.h:418
virtual void Streamer(TBuffer &)
Stream a string object.
Definition TString.cxx:1390
TString & Append(const char *cs)
Definition TString.h:576
void Form(const char *fmt,...)
Formats a string using a printf style format descriptor.
Definition TString.cxx:2334
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:636
Ssiz_t Index(const char *pat, Ssiz_t i=0, ECaseCompare cmp=kExact) const
Definition TString.h:651
Double_t GetTimeOffset() const
Definition TStyle.h:261
Int_t GetOptLogy() const
Definition TStyle.h:240
Int_t GetOptStat() const
Definition TStyle.h:237
Int_t GetOptTitle() const
Definition TStyle.h:238
void SetCanvasBorderSize(Width_t size=1)
Definition TStyle.h:329
Float_t GetScreenFactor() const
Definition TStyle.h:248
Int_t GetPadTickX() const
Definition TStyle.h:209
Bool_t IsReading() const
Definition TStyle.h:283
void SetCanvasColor(Color_t color=19)
Definition TStyle.h:328
Float_t GetPadRightMargin() const
Definition TStyle.h:206
void SetCanvasBorderMode(Int_t mode=1)
Definition TStyle.h:330
Int_t GetCanvasDefH() const
Definition TStyle.h:184
Int_t GetCanvasDefX() const
Definition TStyle.h:186
Bool_t GetPadGridY() const
Definition TStyle.h:208
Float_t GetPadLeftMargin() const
Definition TStyle.h:205
void SetPalette(Int_t ncolors=kBird, Int_t *colors=nullptr, Float_t alpha=1.)
See TColor::SetPalette.
Definition TStyle.cxx:1781
Bool_t GetCanvasPreferGL() const
Definition TStyle.h:180
Int_t GetCanvasDefY() const
Definition TStyle.h:187
Bool_t GetPadGridX() const
Definition TStyle.h:207
Int_t GetPadTickY() const
Definition TStyle.h:210
Color_t GetCanvasColor() const
Definition TStyle.h:181
Float_t GetPadBottomMargin() const
Definition TStyle.h:203
Int_t GetCanvasDefW() const
Definition TStyle.h:185
Int_t GetOptLogx() const
Definition TStyle.h:239
Int_t GetCanvasBorderMode() const
Definition TStyle.h:183
Width_t GetCanvasBorderSize() const
Definition TStyle.h:182
Int_t GetOptFit() const
Definition TStyle.h:236
Int_t GetOptLogz() const
Definition TStyle.h:241
Float_t GetPadTopMargin() const
Definition TStyle.h:204
static TClass * Class()
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.
small helper class to store/restore gPad context in TPad methods
Definition TVirtualPad.h:61
auto GetSaved() const
Definition TVirtualPad.h:69
TVirtualPad is an abstract base class for the Pad and Canvas classes.
Definition TVirtualPad.h:51
static TClass * Class()
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:691
Definition first.py:1
th1 Draw()