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 "TSystem.h"
36#include "TObjArray.h"
37#include "TVirtualPadEditor.h"
38#include "TVirtualViewer3D.h"
39#include "TPadPainter.h"
40#include "TPadPainterPS.h"
41#include "TVirtualGL.h"
42#include "TVirtualPS.h"
43#include "TVirtualX.h"
44#include "TAxis.h"
45#include "TH1.h"
46#include "TGraph.h"
47#include "TMath.h"
48#include "TView.h"
49#include "strlcpy.h"
50#include "snprintf.h"
51
52#include "TVirtualMutex.h"
53
58
59//*-*x16 macros/layout_canvas
60
62
64
65
66TString GetNewCanvasName(const char *arg = nullptr)
67{
68 if (arg && *arg)
69 return arg;
70
71 const char *defcanvas = gROOT->GetDefCanvasName();
73
74 auto lc = (TList*)gROOT->GetListOfCanvases();
75 Int_t n = lc->GetSize() + 1;
76
77 while(lc->FindObject(cdef.Data()))
78 cdef.Form("%s_n%d", defcanvas, n++);
79
80 return cdef;
81}
82
83
84/** \class TCanvas
85\ingroup gpad
86
87The Canvas class.
88
89A Canvas is an area mapped to a window directly under the control of the display
90manager. A ROOT session may have several canvases open at any given time.
91
92A Canvas may be subdivided into independent graphical areas: the __Pads__.
93A canvas has a default pad which has the name of the canvas itself.
94An example of a Canvas layout is sketched in the picture below.
95
96\image html gpad_canvas.png
97
98This canvas contains two pads named P1 and P2. Both Canvas, P1 and P2 can be
99moved, grown, shrunk using the normal rules of the Display manager.
100
101Once objects have been drawn in a canvas, they can be edited/moved by pointing
102directly to them. The cursor shape is changed to suggest the type of action that
103one can do on this object. Clicking with the right mouse button on an object
104pops-up a contextmenu with a complete list of actions possible on this object.
105
106A graphical editor may be started from the canvas "View" menu under the menu
107entry "Toolbar".
108
109An interactive HELP is available by clicking on the HELP button at the top right
110of the canvas. It gives a short explanation about the canvas' menus.
111
112A canvas may be automatically divided into pads via `TPad::Divide`.
113
114At creation time, no matter if in interactive or batch mode, the constructor
115defines the size of the canvas window (including the size of the window
116manager's decoration). To define precisely the graphics area size of a canvas in
117the interactive mode, the following four lines of code should be used:
118~~~ {.cpp}
119 {
120 Double_t w = 600;
121 Double_t h = 600;
122 auto c = new TCanvas("c", "c", w, h);
123 c->SetWindowSize(w + (w - c->GetWw()), h + (h - c->GetWh()));
124 }
125~~~
126and in the batch mode simply do:
127~~~ {.cpp}
128 c->SetCanvasSize(w,h);
129~~~
130
131If the canvas size exceeds the window size, scroll bars will be added to the canvas
132This allows to display very large canvases (even bigger than the screen size). The
133Following example shows how to proceed.
134~~~ {.cpp}
135 {
136 auto c = new TCanvas("c","c");
137 c->SetCanvasSize(1500, 1500);
138 c->SetWindowSize(500, 500);
139 }
140~~~
141*/
142
143////////////////////////////////////////////////////////////////////////////////
144/// Canvas default constructor.
145
146TCanvas::TCanvas(Bool_t build) : TPad(), fDoubleBuffer(0)
147{
148 fPainter = nullptr;
149 fWindowTopX = 0;
150 fWindowTopY = 0;
151 fWindowWidth = 0;
152 fWindowHeight = 0;
153 fCw = 0;
154 fCh = 0;
155 fXsizeUser = 0;
156 fYsizeUser = 0;
159 fHighLightColor = gEnv->GetValue("Canvas.HighLightColor", kRed);
160 fEvent = -1;
161 fEventX = -1;
162 fEventY = -1;
163 fSelectedX = 0;
164 fSelectedY = 0;
166 fDrawn = kFALSE;
168 fSelected = nullptr;
169 fClickSelected = nullptr;
170 fSelectedPad = nullptr;
171 fClickSelectedPad = nullptr;
172 fPadSave = nullptr;
173 fCanvasImp = nullptr;
174 fContextMenu = nullptr;
175
177
178 if (!build || TClass::IsCallingNew() != TClass::kRealNew) {
179 Constructor();
180 } else {
181 auto cdef = GetNewCanvasName();
182
183 Constructor(cdef.Data(), cdef.Data(), 1);
184 }
185}
186
187////////////////////////////////////////////////////////////////////////////////
188/// Canvas default constructor
189
191{
192 if (gThreadXAR) {
193 void *arr[2];
194 arr[1] = this;
195 if ((*gThreadXAR)("CANV", 2, arr, nullptr)) return;
196 }
197
198 fCanvas = nullptr;
199 fCanvasID = -1;
200 fCanvasImp = nullptr;
201 fBatch = kTRUE;
203
204 fContextMenu = nullptr;
205 fSelected = nullptr;
206 fClickSelected = nullptr;
207 fSelectedPad = nullptr;
208 fClickSelectedPad = nullptr;
209 fPadSave = nullptr;
213}
214
215////////////////////////////////////////////////////////////////////////////////
216/// Create an embedded canvas, i.e. a canvas that is in a TGCanvas widget
217/// which is placed in a TGFrame. This ctor is only called via the
218/// TRootEmbeddedCanvas class.
219///
220/// If "name" starts with "gl" the canvas is ready to receive GL output.
221
222TCanvas::TCanvas(const char *name, Int_t ww, Int_t wh, Int_t winid) : TPad(), fDoubleBuffer(0)
223{
224 fCanvasImp = nullptr;
225 fPainter = nullptr;
226 Init();
227
229 fWindowTopX = 0;
230 fWindowTopY = 0;
231 fWindowWidth = ww;
233 fCw = ww + 4;
234 fCh = wh +28;
235 fBatch = kFALSE;
237
238 //This is a very special ctor. A window exists already!
239 //Can create painter now.
241
242 if (fUseGL) {
243 fGLDevice = gGLManager->CreateGLContext(winid);
244 if (fGLDevice == -1)
245 fUseGL = kFALSE;
246 }
247
249 if (!fCanvasImp) return;
250
252 fName = GetNewCanvasName(name); // avoid Modified() signal from SetName
253 Build();
254}
255
256////////////////////////////////////////////////////////////////////////////////
257/// Create a new canvas with a predefined size form.
258/// If form < 0 the menubar is not shown.
259///
260/// - form = 1 700x500 at 10,10 (set by TStyle::SetCanvasDefH,W,X,Y)
261/// - form = 2 500x500 at 20,20
262/// - form = 3 500x500 at 30,30
263/// - form = 4 500x500 at 40,40
264/// - form = 5 500x500 at 50,50
265///
266/// If "name" starts with "gl" the canvas is ready to receive GL output.
267
268TCanvas::TCanvas(const char *name, const char *title, Int_t form) : TPad(), fDoubleBuffer(0)
269{
270 fPainter = nullptr;
272
273 Constructor(name, title, form);
274}
275
276////////////////////////////////////////////////////////////////////////////////
277/// Create a new canvas with a predefined size form.
278/// If form < 0 the menubar is not shown.
279///
280/// - form = 1 700x500 at 10,10 (set by TStyle::SetCanvasDefH,W,X,Y)
281/// - form = 2 500x500 at 20,20
282/// - form = 3 500x500 at 30,30
283/// - form = 4 500x500 at 40,40
284/// - form = 5 500x500 at 50,50
285
286void TCanvas::Constructor(const char *name, const char *title, Int_t form)
287{
288 if (gThreadXAR) {
289 void *arr[6];
290 static Int_t ww = 500;
291 static Int_t wh = 500;
292 arr[1] = this; arr[2] = (void*)name; arr[3] = (void*)title; arr[4] =&ww; arr[5] = &wh;
293 if ((*gThreadXAR)("CANV", 6, arr, nullptr)) return;
294 }
295
296 Init();
297 SetBit(kMenuBar,true);
298 if (form < 0) {
299 form = -form;
300 SetBit(kMenuBar,false);
301 }
302
303 fCanvas = this;
304
305 fCanvasID = -1;
306 TCanvas *old = (TCanvas*)gROOT->GetListOfCanvases()->FindObject(name);
307 if (old && old->IsOnHeap()) {
308 Warning("Constructor","Deleting canvas with same name: %s",name);
309 delete old;
310 }
311 if (gROOT->IsBatch()) { //We are in Batch mode
313 if (form == 1) {
316 } else {
317 fWindowWidth = 500;
318 fWindowHeight = 500;
319 }
323 if (!fCanvasImp) return;
324 fBatch = kTRUE;
325 } else { //normal mode with a screen window
327 if (form < 1 || form > 20) form = 1;
328 auto factory = gROOT->IsWebDisplay() ? gBatchGuiFactory : gGuiFactory;
329 Int_t ux, uy, cw, ch;
330 if (form == 1) {
331 cw = gStyle->GetCanvasDefW();
332 ch = gStyle->GetCanvasDefH();
335 } else {
336 cw = ch = 500;
337 ux = uy = form * 10;
338 }
339
340 fCanvasImp = factory->CreateCanvasImp(this, name, Int_t(cx*ux), Int_t(cx*uy), UInt_t(cx*cw), UInt_t(cx*ch));
341 if (!fCanvasImp) return;
342
343 if (!gROOT->IsBatch() && fCanvasID == -1)
345
347 fBatch = kFALSE;
348 }
349
351
352 fName = GetNewCanvasName(name); // avoid Modified() signal from SetName
353 SetTitle(title); // requires fCanvasImp set
354 Build();
355
356 // Popup canvas
357 fCanvasImp->Show();
358}
359
360////////////////////////////////////////////////////////////////////////////////
361/// Create a new canvas at a random position.
362///
363/// \param[in] name canvas name
364/// \param[in] title canvas title
365/// \param[in] ww is the window size in pixels along X
366/// (if ww < 0 the menubar is not shown)
367/// \param[in] wh is the window size in pixels along Y
368///
369/// If "name" starts with "gl" the canvas is ready to receive GL output.
370
371TCanvas::TCanvas(const char *name, const char *title, Int_t ww, Int_t wh) : TPad(), fDoubleBuffer(0)
372{
373 fPainter = nullptr;
375
376 Constructor(name, title, ww, wh);
377}
378
379////////////////////////////////////////////////////////////////////////////////
380/// Create a new canvas at a random position.
381///
382/// \param[in] name canvas name
383/// \param[in] title canvas title
384/// \param[in] ww is the window size in pixels along X
385/// (if ww < 0 the menubar is not shown)
386/// \param[in] wh is the window size in pixels along Y
387
388void TCanvas::Constructor(const char *name, const char *title, Int_t ww, Int_t wh)
389{
390 if (gThreadXAR) {
391 void *arr[6];
392 arr[1] = this; arr[2] = (void*)name; arr[3] = (void*)title; arr[4] =&ww; arr[5] = &wh;
393 if ((*gThreadXAR)("CANV", 6, arr, nullptr)) return;
394 }
395
396 Init();
397 SetBit(kMenuBar,true);
398 if (ww < 0) {
399 ww = -ww;
400 SetBit(kMenuBar,false);
401 }
402 if (wh <= 0) {
403 Error("Constructor", "Invalid canvas height: %d",wh);
404 return;
405 }
406 fCw = ww;
407 fCh = wh;
408 fCanvasID = -1;
409 TCanvas *old = (TCanvas*)gROOT->GetListOfCanvases()->FindObject(name);
410 if (old && old->IsOnHeap()) {
411 Warning("Constructor","Deleting canvas with same name: %s",name);
412 delete old;
413 }
414 if (gROOT->IsBatch()) { //We are in Batch mode
416 fWindowWidth = ww;
418 fCw = ww;
419 fCh = wh;
421 if (!fCanvasImp) return;
422 fBatch = kTRUE;
423 } else {
425 auto factory = gROOT->IsWebDisplay() ? gBatchGuiFactory : gGuiFactory;
426 fCanvasImp = factory->CreateCanvasImp(this, name, UInt_t(cx*ww), UInt_t(cx*wh));
427 if (!fCanvasImp) return;
428
429 if (!gROOT->IsBatch() && fCanvasID == -1)
431
433 fBatch = kFALSE;
434 }
435
437
438 fName = GetNewCanvasName(name); // avoid Modified() signal from SetName
439 SetTitle(title); // requires fCanvasImp set
440 Build();
441
442 // Popup canvas
443 fCanvasImp->Show();
444}
445
446////////////////////////////////////////////////////////////////////////////////
447/// Create a new canvas.
448///
449/// \param[in] name canvas name
450/// \param[in] title canvas title
451/// \param[in] wtopx,wtopy are the pixel coordinates of the top left corner of
452/// the canvas (if wtopx < 0) the menubar is not shown)
453/// \param[in] ww is the window size in pixels along X
454/// \param[in] wh is the window size in pixels along Y
455///
456/// If "name" starts with "gl" the canvas is ready to receive GL output.
457
458TCanvas::TCanvas(const char *name, const char *title, Int_t wtopx, Int_t wtopy, Int_t ww, Int_t wh)
459 : TPad(), fDoubleBuffer(0)
460{
461 fPainter = nullptr;
463
464 Constructor(name, title, wtopx, wtopy, ww, wh);
465}
466
467////////////////////////////////////////////////////////////////////////////////
468/// Create a new canvas.
469///
470/// \param[in] name canvas name
471/// \param[in] title canvas title
472/// \param[in] wtopx,wtopy are the pixel coordinates of the top left corner of
473/// the canvas (if wtopx < 0) the menubar is not shown)
474/// \param[in] ww is the window size in pixels along X
475/// \param[in] wh is the window size in pixels along Y
476
477void TCanvas::Constructor(const char *name, const char *title, Int_t wtopx,
478 Int_t wtopy, Int_t ww, Int_t wh)
479{
480 if (gThreadXAR) {
481 void *arr[8];
482 arr[1] = this; arr[2] = (void*)name; arr[3] = (void*)title;
483 arr[4] = &wtopx; arr[5] = &wtopy; arr[6] = &ww; arr[7] = &wh;
484 if ((*gThreadXAR)("CANV", 8, arr, nullptr)) return;
485 }
486
487 Init();
488 SetBit(kMenuBar,true);
489 if (wtopx < 0) {
490 wtopx = -wtopx;
491 SetBit(kMenuBar,false);
492 }
493 fCw = ww;
494 fCh = wh;
495 fCanvasID = -1;
496 TCanvas *old = (TCanvas*)gROOT->GetListOfCanvases()->FindObject(name);
497 if (old && old->IsOnHeap()) {
498 Warning("Constructor","Deleting canvas with same name: %s",name);
499 delete old;
500 }
501 if (gROOT->IsBatch()) { //We are in Batch mode
503 fWindowWidth = ww;
505 fCw = ww;
506 fCh = wh;
508 if (!fCanvasImp) return;
509 fBatch = kTRUE;
510 } else { //normal mode with a screen window
512 auto factory = gROOT->IsWebDisplay() ? gBatchGuiFactory : gGuiFactory;
513 fCanvasImp = factory->CreateCanvasImp(this, name, Int_t(cx*wtopx), Int_t(cx*wtopy), UInt_t(cx*ww), UInt_t(cx*wh));
514 if (!fCanvasImp) return;
515
516 if (!gROOT->IsBatch() && fCanvasID == -1)
518
520 fBatch = kFALSE;
521 }
522
524
525 fName = GetNewCanvasName(name); // avoid Modified() signal from SetName
526 SetTitle(title); // requires fCanvasImp set
527 Build();
528
529 // Popup canvas
530 fCanvasImp->Show();
531}
532
533////////////////////////////////////////////////////////////////////////////////
534/// Initialize the TCanvas members. Called by all constructors.
535
537{
538 // Make sure the application environment exists. It is need for graphics
539 // (colors are initialized in the TApplication ctor).
540 if (!gApplication)
542
543 // Load and initialize graphics libraries if
544 // TApplication::NeedGraphicsLibs() has been called by a
545 // library static initializer.
546 if (gApplication)
547 gApplication->InitializeGraphics(gROOT->IsWebDisplay());
548
549 // Get some default from .rootrc. Used in fCanvasImp->InitWindow().
550 fHighLightColor = gEnv->GetValue("Canvas.HighLightColor", kRed);
551 SetBit(kMoveOpaque, gEnv->GetValue("Canvas.MoveOpaque", 0));
552 SetBit(kResizeOpaque, gEnv->GetValue("Canvas.ResizeOpaque", 0));
553 if (gEnv->GetValue("Canvas.ShowEventStatus", kFALSE)) SetBit(kShowEventStatus);
554 if (gEnv->GetValue("Canvas.ShowToolTips", kFALSE)) SetBit(kShowToolTips);
555 if (gEnv->GetValue("Canvas.ShowToolBar", kFALSE)) SetBit(kShowToolBar);
556 if (gEnv->GetValue("Canvas.ShowEditor", kFALSE)) SetBit(kShowEditor);
557 if (gEnv->GetValue("Canvas.AutoExec", kTRUE)) SetBit(kAutoExec);
558
559 // Fill canvas ROOT data structure
560 fXsizeUser = 0;
561 fYsizeUser = 0;
564
565 fDISPLAY = "$DISPLAY";
568 fSelected = nullptr;
569 fClickSelected = nullptr;
570 fSelectedX = 0;
571 fSelectedY = 0;
572 fSelectedPad = nullptr;
573 fClickSelectedPad= nullptr;
574 fPadSave = nullptr;
575 fEvent = -1;
576 fEventX = -1;
577 fEventY = -1;
578 fContextMenu = nullptr;
579 fDrawn = kFALSE;
581}
582
583////////////////////////////////////////////////////////////////////////////////
584/// Build a canvas. Called by all constructors.
585
587{
588 // Get window identifier
589 if (fCanvasID == -1 && fCanvasImp)
591 if (fCanvasID == -1) return;
592
593 if (fCw !=0 && fCh !=0) {
596 }
597
598 // Set Pad parameters
599 gPad = this;
600 fCanvas = this;
601 fMother = (TPad*)gPad;
602
603 if (IsBatch()) {
604 // Make sure that batch interactive canvas sizes are the same
605 fCw -= 4;
606 fCh -= 28;
607 } else if (IsWeb()) {
608 // mark canvas as batch - avoid gVirtualX in many places
610 } else {
611 //normal mode with a screen window
612 // Set default physical canvas attributes
614 fPainter->SetAttFill({1, 1001}); //Set color index for fill area
615 fPainter->SetAttLine({1, 1, 1}); //Set color index for lines
616 fPainter->SetAttMarker({1, 1, 1}); //Set color index for markers
617 fPainter->SetAttText({22, 0., 1, 42, 12}); //Set color index for text
618 // Clear workstation
620
621 // Set Double Buffer on by default
623
624 // Get effective window parameters (with borders and menubar)
627
628 // Get effective canvas parameters without borders
630
631 fContextMenu = new TContextMenu("ContextMenu");
632 }
633
634 gROOT->GetListOfCanvases()->Add(this);
635
636 if (!fPrimitives) {
637 fPrimitives = new TList;
639 SetFillStyle(1001);
651 fBorderMode=gStyle->GetCanvasBorderMode(); // do not call SetBorderMode (function redefined in TCanvas)
652 SetPad(0, 0, 1, 1);
653 Range(0, 0, 1, 1); //pad range is set by default to [0,1] in x and y
654
656 if (vpp) vpp->SelectDrawable(fPixmapID);//gVirtualX->SelectPixmap(fPixmapID); //pixmap must be selected
657 PaintBorder(GetFillColor(), kTRUE); //paint background
658 }
659
660 // transient canvases have typically no menubar and should not get
661 // by default the event status bar (if set by default)
662 if (TestBit(kMenuBar) && fCanvasImp) {
664 // ... and toolbar + editor
668 }
669}
670
671////////////////////////////////////////////////////////////////////////////////
672/// Canvas destructor
673
675{
676 Destructor();
677}
678
679////////////////////////////////////////////////////////////////////////////////
680/// Browse.
681
683{
684 Draw();
685 cd();
687}
688
689////////////////////////////////////////////////////////////////////////////////
690/// Actual canvas destructor.
691
693{
694 if (gThreadXAR) {
695 void *arr[2];
696 arr[1] = this;
697 if ((*gThreadXAR)("CDEL", 2, arr, nullptr)) return;
698 }
699
700 if (ROOT::Detail::HasBeenDeleted(this)) return;
701
703 if (!gPad) return;
704
705 Close();
706
707 //If not yet (batch mode?).
709}
710
711////////////////////////////////////////////////////////////////////////////////
712/// Set current canvas & pad. Returns the new current pad,
713/// or 0 in case of failure.
714/// See TPad::cd() for an explanation of the parameter.
715
717{
718 if (fCanvasID == -1) return nullptr;
719
721
722 // in case doublebuffer is off, draw directly onto display window
723 if (!IsBatch() && !IsWeb() && !fDoubleBuffer && fPainter)
724 fPainter->SelectDrawable(fCanvasID);//Ok, does not matter for glpad.
725
726 return gPad;
727}
728
729////////////////////////////////////////////////////////////////////////////////
730/// Remove all primitives from the canvas.
731/// If option "D" is specified, direct sub-pads are cleared but not deleted.
732/// This option is not recursive, i.e. pads in direct sub-pads are deleted.
733
735{
736 if (fCanvasID == -1) return;
737
739
740 TString opt = option;
741 opt.ToLower();
742 if (opt.Contains("d")) {
743 // clear subpads, but do not delete pads in case the canvas
744 // has been divided (note: option "D" is propagated so could cause
745 // conflicts for primitives using option "D" for something else)
746 if (fPrimitives) {
747 TIter next(fPrimitives);
748 TObject *obj;
749 while ((obj=next())) {
750 obj->Clear(option);
751 }
752 }
753 } else {
754 //default, clear everything in the canvas. Subpads are deleted
755 TPad::Clear(option); //Remove primitives from pad
756 }
757
758 fSelected = nullptr;
759 fClickSelected = nullptr;
760 fSelectedPad = nullptr;
761 fClickSelectedPad = nullptr;
762}
763
764////////////////////////////////////////////////////////////////////////////////
765/// Emit pad Cleared signal.
766
768{
769 Emit("Cleared(TVirtualPad*)", (Longptr_t)pad);
770}
771
772////////////////////////////////////////////////////////////////////////////////
773/// Emit Closed signal.
774
776{
777 Emit("Closed()");
778}
779
780////////////////////////////////////////////////////////////////////////////////
781/// Close canvas.
782///
783/// Delete window/pads data structure
784
786{
787 auto padsave = gPad;
788 TCanvas *cansave = padsave ? padsave->GetCanvas() : nullptr;
789
790 if (fCanvasID != -1) {
791
792 if (!gROOT->IsLineProcessing() && !gVirtualX->IsCmdThread()) {
793 gInterpreter->Execute(this, IsA(), "Close", option);
794 return;
795 }
796
798
800
801 cd();
803
804 if (!IsBatch() && !IsWeb()) {
805 //select current canvas
806 if (fPainter)
808
810
811 if (fCanvasImp)
812 fCanvasImp->Close();
813 }
814 fCanvasID = -1;
815 fBatch = kTRUE;
816
817 gROOT->GetListOfCanvases()->Remove(this);
818
819 // Close actual window on screen
821 }
822
823 if (cansave == this) {
824 gPad = (TCanvas *) gROOT->GetListOfCanvases()->First();
825 } else {
826 gPad = padsave;
827 }
828
829 Closed();
830}
831
832////////////////////////////////////////////////////////////////////////////////
833/// Copy the canvas pixmap of the pad to the canvas.
834
836{
837 if (!IsBatch()) {
840 }
841}
842
843////////////////////////////////////////////////////////////////////////////////
844/// Draw a canvas.
845/// If a canvas with the name is already on the screen, the canvas is repainted.
846/// This function is useful when a canvas object has been saved in a Root file.
847/// One can then do:
848/// ~~~ {.cpp}
849/// Root > TFile::Open("file.root");
850/// Root > canvas->Draw();
851/// ~~~
852
854{
855 // Load and initialize graphics libraries if
856 // TApplication::NeedGraphicsLibs() has been called by a
857 // library static initializer.
858 if (gApplication)
859 gApplication->InitializeGraphics(gROOT->IsWebDisplay());
860
861 fDrawn = kTRUE;
862
863 TCanvas *old = (TCanvas*)gROOT->GetListOfCanvases()->FindObject(GetName());
864 if (old == this) {
865 if (IsWeb()) {
866 Modified();
867 UpdateAsync();
868 } else {
869 Paint();
870 }
871 return;
872 }
873 if (old) { gROOT->GetListOfCanvases()->Remove(old); delete old;}
874
875 if (fWindowWidth == 0) {
876 if (fCw !=0) fWindowWidth = fCw+4;
877 else fWindowWidth = 800;
878 }
879 if (fWindowHeight == 0) {
880 if (fCh !=0) fWindowHeight = fCh+28;
881 else fWindowHeight = 600;
882 }
883 if (gROOT->IsBatch()) { //We are in Batch mode
885 if (!fCanvasImp) return;
886 fBatch = kTRUE;
887
888 } else { //normal mode with a screen window
889 auto factory = gROOT->IsWebDisplay() ? gBatchGuiFactory : gGuiFactory;
890 fCanvasImp = factory->CreateCanvasImp(this, GetName(), fWindowTopX, fWindowTopY,
892 if (!fCanvasImp) return;
894 }
895 Build();
896 ResizePad();
898 fCanvasImp->Show();
899 Modified();
900}
901
902////////////////////////////////////////////////////////////////////////////////
903/// Draw a clone of this canvas
904/// A new canvas is created that is a clone of this canvas
905
907{
909 newCanvas->SetName();
910
911 newCanvas->Draw(option);
912 newCanvas->Update();
913 return newCanvas;
914}
915
916////////////////////////////////////////////////////////////////////////////////
917/// Draw a clone of this canvas into the current pad
918/// In an interactive session, select the destination/current pad
919/// with the middle mouse button, then point to the canvas area to select
920/// the canvas context menu item DrawClonePad.
921/// Note that the original canvas may have subpads.
922
924{
925 auto padsav = gPad;
926 auto selpad = gROOT->GetSelectedPad();
927 auto pad = padsav;
928 if (pad == this) pad = selpad;
929 if (!padsav || !pad || pad == this) {
931 newCanvas->SetWindowSize(GetWindowWidth(),GetWindowHeight());
932 return newCanvas;
933 }
934 if (fCanvasID == -1) {
935 auto factory = gROOT->IsWebDisplay() ? gBatchGuiFactory : gGuiFactory;
936 fCanvasImp = factory->CreateCanvasImp(this, GetName(), fWindowTopX, fWindowTopY,
938 if (!fCanvasImp) return nullptr;
941 }
942 this->cd();
943 //copy pad attributes
944 pad->Range(fX1,fY1,fX2,fY2);
945 pad->SetTickx(GetTickx());
946 pad->SetTicky(GetTicky());
947 pad->SetGridx(GetGridx());
948 pad->SetGridy(GetGridy());
949 pad->SetLogx(GetLogx());
950 pad->SetLogy(GetLogy());
951 pad->SetLogz(GetLogz());
952 pad->SetBorderSize(GetBorderSize());
953 pad->SetBorderMode(GetBorderMode());
957
958 //copy primitives
960 while (auto obj = next()) {
961 pad->cd();
962 pad->Add(obj->Clone(), next.GetOption(), kFALSE); // do not issue modified for each object
963 }
964 pad->ResizePad();
965 pad->Modified();
966 pad->Update();
967 if (padsav) padsav->cd();
968 return nullptr;
969}
970
971////////////////////////////////////////////////////////////////////////////////
972/// Report name and title of primitive below the cursor.
973///
974/// This function is called when the option "Event Status"
975/// in the canvas menu "Options" is selected.
976
977void TCanvas::DrawEventStatus(Int_t event, Int_t px, Int_t py, TObject *selected)
978{
979 const Int_t kTMAX=256;
980 static char atext[kTMAX];
981
982 if (!TestBit(kShowEventStatus) || !selected) return;
983
984 if (!fCanvasImp) return; //this may happen when closing a TAttCanvas
985
987
988 fCanvasImp->SetStatusText(selected->GetTitle(),0);
989 fCanvasImp->SetStatusText(selected->GetName(),1);
990 if (event == kKeyPress)
991 snprintf(atext, kTMAX, "%c", (char) px);
992 else
993 snprintf(atext, kTMAX, "%d,%d", px, py);
995
996 // Show date/time if TimeDisplay is selected
997 TAxis *xaxis = nullptr;
998 if ( selected->InheritsFrom("TH1") )
999 xaxis = ((TH1*)selected)->GetXaxis();
1000 else if ( selected->InheritsFrom("TGraph") )
1001 xaxis = ((TGraph*)selected)->GetXaxis();
1002 else if ( selected->InheritsFrom("TAxis") )
1003 xaxis = (TAxis*)selected;
1004 if ( xaxis != nullptr && xaxis->GetTimeDisplay()) {
1005 TString objinfo = selected->GetObjectInfo(px,py);
1006 // check if user has overwritten GetObjectInfo and altered
1007 // the default text from TObject::GetObjectInfo "x=.. y=.."
1008 if (objinfo.Contains("x=") && objinfo.Contains("y=") ) {
1009 UInt_t toff = 0;
1010 TString time_format(xaxis->GetTimeFormat());
1011 // TimeFormat may contain offset: %F2000-01-01 00:00:00
1012 Int_t idF = time_format.Index("%F");
1013 if (idF>=0) {
1014 Int_t lnF = time_format.Length();
1015 // minimal check for correct format
1016 if (lnF - idF == 21) {
1019 toff = dtoff.Convert();
1020 }
1021 } else {
1023 }
1024 TDatime dt((UInt_t)gPad->AbsPixeltoX(px) + toff);
1025 snprintf(atext, kTMAX, "%s, y=%g",
1026 dt.AsSQLString(),gPad->AbsPixeltoY(py));
1028 return;
1029 }
1030 }
1031 // default
1032 fCanvasImp->SetStatusText(selected->GetObjectInfo(px,py),3);
1033}
1034
1035////////////////////////////////////////////////////////////////////////////////
1036/// Get editor bar.
1037
1042
1043////////////////////////////////////////////////////////////////////////////////
1044/// Embedded a canvas into a TRootEmbeddedCanvas. This method is only called
1045/// via TRootEmbeddedCanvas::AdoptCanvas.
1046
1048{
1049 // If fCanvasImp already exists, no need to go further.
1050 if(fCanvasImp) return;
1051
1052 fCanvasID = winid;
1053 fWindowTopX = 0;
1054 fWindowTopY = 0;
1055 fWindowWidth = ww;
1056 fWindowHeight = wh;
1057 fCw = ww;
1058 fCh = wh;
1059 fBatch = kFALSE;
1060 fUpdating = kFALSE;
1061
1063 if (!fCanvasImp) return;
1064 Build();
1065 Resize();
1066}
1067
1068////////////////////////////////////////////////////////////////////////////////
1069/// Generate kMouseEnter and kMouseLeave events depending on the previously
1070/// selected object and the currently selected object. Does nothing if the
1071/// selected object does not change.
1072
1074{
1075 if (prevSelObj == fSelected) return;
1076
1079
1080 if (prevSelObj) {
1081 gPad = prevSelPad;
1082 prevSelObj->ExecuteEvent(kMouseLeave, fEventX, fEventY);
1084 RunAutoExec();
1086 }
1087
1089
1090 if (fSelected) {
1093 RunAutoExec();
1095 }
1096
1097 fEvent = sevent;
1098}
1099
1100////////////////////////////////////////////////////////////////////////////////
1101/// Execute action corresponding to one event.
1102///
1103/// This member function must be implemented to realize the action
1104/// corresponding to the mouse click on the object in the canvas
1105///
1106/// Only handle mouse motion events in TCanvas, all other events are
1107/// ignored for the time being
1108
1110{
1111 if (gROOT->GetEditorMode()) {
1112 TPad::ExecuteEvent(event,px,py);
1113 return;
1114 }
1115
1116 switch (event) {
1117
1118 case kMouseMotion:
1120 break;
1121 }
1122}
1123
1124////////////////////////////////////////////////////////////////////////////////
1125/// Turn rubberband feedback mode on or off.
1126
1128{
1129 if (IsWeb() || (fCanvasID == -1))
1130 return;
1131
1132 SetDoubleBuffer(set ? 0 : 1); // switch double buffer
1133
1134 if (fPainter)
1136}
1137
1138////////////////////////////////////////////////////////////////////////////////
1139/// Flush canvas buffers.
1140
1142{
1143 if ((fCanvasID == -1) || IsWeb()) return;
1144
1145 TContext ctxt(this, kTRUE);
1146 if (!IsBatch()) {
1147 if (!UseGL() || fGLDevice == -1) {
1149 gPad = ctxt.GetSaved(); //don't do cd() because than also the pixmap is changed
1150 CopyPixmaps();
1152 } else {
1154 gVirtualPS = nullptr;
1155 gGLManager->MakeCurrent(fGLDevice);
1157 Paint();
1158 if (ctxt.GetSaved() && ctxt.GetSaved()->GetCanvas() == this) {
1159 ctxt.GetSaved()->cd();
1160 ctxt.GetSaved()->HighLight(ctxt.GetSaved()->GetHighLightColor());
1161 //cd();
1162 }
1164 gGLManager->Flush(fGLDevice);
1165 gVirtualPS = tvps;
1166 }
1167 }
1168}
1169
1170////////////////////////////////////////////////////////////////////////////////
1171/// Force canvas update
1172
1177
1178////////////////////////////////////////////////////////////////////////////////
1179/// Force a copy of current style for all objects in canvas.
1180
1182{
1183 if (!gROOT->IsLineProcessing() && !gVirtualX->IsCmdThread()) {
1184 gInterpreter->Execute(this, IsA(), "UseCurrentStyle", "");
1185 return;
1186 }
1187
1189
1191
1192 if (gStyle->IsReading()) {
1196 } else {
1200 }
1201}
1202
1203////////////////////////////////////////////////////////////////////////////////
1204/// Returns current top x position of window on screen.
1205
1213
1214////////////////////////////////////////////////////////////////////////////////
1215/// Returns current top y position of window on screen.
1216
1224
1225////////////////////////////////////////////////////////////////////////////////
1226/// Handle Input Events.
1227///
1228/// Handle input events, like button up/down in current canvas.
1229
1231{
1232 TPad *pad;
1235
1236 fPadSave = (TPad*)gPad;
1237 cd(); // make sure this canvas is the current canvas
1238
1239 fEvent = event;
1240 fEventX = px;
1241 fEventY = py;
1242
1243 switch (event) {
1244
1245 case kMouseMotion:
1246 // highlight object tracked over
1247 pad = Pick(px, py, prevSelObj);
1248 if (!pad) return;
1249
1251
1252 gPad = pad; // don't use cd() we will use the current
1253 // canvas via the GetCanvas member and not via
1254 // gPad->GetCanvas
1255
1256 if (fSelected) {
1257 fSelected->ExecuteEvent(event, px, py);
1258 RunAutoExec();
1259 }
1260
1261 break;
1262
1263 case kMouseEnter:
1264 // mouse enters canvas
1266 break;
1267
1268 case kMouseLeave:
1269 // mouse leaves canvas
1270 {
1271 // force popdown of tooltips
1274 fSelected = nullptr;
1275 fSelectedPad = nullptr;
1277 fSelected = sobj;
1280 }
1281 break;
1282
1283 case kButton1Double:
1284 // triggered on the second button down within 350ms and within
1285 // 3x3 pixels of the first button down, button up finishes action
1286
1287 case kButton1Down:
1288 // find pad in which input occurred
1289 pad = Pick(px, py, prevSelObj);
1290 if (!pad) return;
1291
1292 gPad = pad; // don't use cd() because we won't draw in pad
1293 // we will only use its coordinate system
1294
1295 if (fSelected) {
1296 FeedbackMode(kTRUE); // to draw in rubberband mode
1297 fSelected->ExecuteEvent(event, px, py);
1298
1299 RunAutoExec();
1300 }
1301
1302 break;
1303
1304 case kArrowKeyPress:
1305 case kArrowKeyRelease:
1306 case kButton1Motion:
1307 case kButton1ShiftMotion: //8 == kButton1Motion + shift modifier
1308 if (fSelected) {
1310
1311 fSelected->ExecuteEvent(event, px, py);
1314 Bool_t resize = kFALSE;
1316 resize = ((TBox*)fSelected)->IsBeingResized();
1318 resize = ((TVirtualPad*)fSelected)->IsBeingResized();
1319
1320 if ((!resize && TestBit(kMoveOpaque)) || (resize && TestBit(kResizeOpaque))) {
1321 gPad = fPadSave;
1322 Update();
1324 }
1325 }
1326
1327 RunAutoExec();
1328 }
1329
1330 break;
1331
1332 case kButton1Up:
1333
1334 if (fSelected) {
1336
1337 fSelected->ExecuteEvent(event, px, py);
1338
1339 RunAutoExec();
1340
1341 if (fPadSave)
1342 gPad = fPadSave;
1343 else {
1344 gPad = this;
1345 fPadSave = this;
1346 }
1347
1348 Update(); // before calling update make sure gPad is reset
1349 }
1350 break;
1351
1352//*-*----------------------------------------------------------------------
1353
1354 case kButton2Down:
1355 // find pad in which input occurred
1356 pad = Pick(px, py, prevSelObj);
1357 if (!pad) return;
1358
1359 gPad = pad; // don't use cd() because we won't draw in pad
1360 // we will only use its coordinate system
1361
1363
1364 if (fSelected) fSelected->Pop(); // pop object to foreground
1365 pad->cd(); // and make its pad the current pad
1366 if (gDebug)
1367 printf("Current Pad: %s / %s\n", pad->GetName(), pad->GetTitle());
1368
1369 // loop over all canvases to make sure that only one pad is highlighted
1370 {
1371 TIter next(gROOT->GetListOfCanvases());
1372 TCanvas *tc;
1373 while ((tc = (TCanvas *)next()))
1374 tc->Update();
1375 }
1376
1377 //if (pad->GetGLDevice() != -1 && fSelected)
1378 // fSelected->ExecuteEvent(event, px, py);
1379
1380 break; // don't want fPadSave->cd() to be executed at the end
1381
1382 case kButton2Motion:
1383 //was empty!
1384 case kButton2Up:
1385 if (fSelected) {
1387
1388 fSelected->ExecuteEvent(event, px, py);
1389 RunAutoExec();
1390 }
1391 break;
1392
1393 case kButton2Double:
1394 break;
1395
1396//*-*----------------------------------------------------------------------
1397
1398 case kButton3Down:
1399 // popup context menu
1400 pad = Pick(px, py, prevSelObj);
1401 if (!pad) return;
1402
1404
1406 !pad->TestBit(kNoContextMenu) && !TestBit(kNoContextMenu))
1407 fContextMenu->Popup(px, py, fSelected, this, pad);
1408
1409 break;
1410
1411 case kButton3Motion:
1412 break;
1413
1414 case kButton3Up:
1416 break;
1417
1418 case kButton3Double:
1419 break;
1420
1421 case kKeyPress:
1422 if (!fSelectedPad || !fSelected) return;
1423 gPad = fSelectedPad; // don't use cd() because we won't draw in pad
1424 // we will only use its coordinate system
1425 fSelected->ExecuteEvent(event, px, py);
1426
1427 RunAutoExec();
1428
1429 break;
1430
1431 case kButton1Shift:
1432 // Try to select
1433 pad = Pick(px, py, prevSelObj);
1434
1435 if (!pad) return;
1436
1438
1439 gPad = pad; // don't use cd() we will use the current
1440 // canvas via the GetCanvas member and not via
1441 // gPad->GetCanvas
1442 if (fSelected) {
1443 fSelected->ExecuteEvent(event, px, py);
1444 RunAutoExec();
1445 }
1446 break;
1447
1448 case kWheelUp:
1449 case kWheelDown:
1450 pad = Pick(px, py, prevSelObj);
1451 if (!pad) return;
1452
1453 gPad = pad;
1454 if (fSelected)
1455 fSelected->ExecuteEvent(event, px, py);
1456 break;
1457
1458 default:
1459 break;
1460 }
1461
1462 if (fPadSave && event != kButton2Down)
1463 fPadSave->cd();
1464
1465 if (event != kMouseLeave) { // signal was already emitted for this event
1466 ProcessedEvent(event, px, py, fSelected); // emit signal
1467 DrawEventStatus(event, px, py, fSelected);
1468 }
1469}
1470
1471////////////////////////////////////////////////////////////////////////////////
1472/// Iconify canvas
1473
1475{
1476 if (fCanvasImp)
1478}
1479
1480////////////////////////////////////////////////////////////////////////////////
1481/// Is folder ?
1482
1484{
1485 return fgIsFolder;
1486}
1487
1488////////////////////////////////////////////////////////////////////////////////
1489/// Is web canvas
1490
1492{
1493 return fCanvasImp ? fCanvasImp->IsWeb() : kFALSE;
1494}
1495
1496////////////////////////////////////////////////////////////////////////////////
1497/// List all pads.
1498
1500{
1502 std::cout <<"Canvas Name=" <<GetName()<<" Title="<<GetTitle()<<" Option="<<option<<std::endl;
1506}
1507
1508////////////////////////////////////////////////////////////////////////////////
1509/// Static function to build a default canvas.
1510
1512{
1513 auto cdef = GetNewCanvasName();
1514
1515 auto c = new TCanvas(cdef.Data(), cdef.Data(), 1);
1516
1517 ::Info("TCanvas::MakeDefCanvas"," created default TCanvas with name %s", cdef.Data());
1518 return c;
1519}
1520
1521////////////////////////////////////////////////////////////////////////////////
1522/// Set option to move objects/pads in a canvas.
1523///
1524/// - set = 1 (default) graphics objects are moved in opaque mode
1525/// - set = 0 only the outline of objects is drawn when moving them
1526///
1527/// The option opaque produces the best effect. It requires however a
1528/// a reasonably fast workstation or response time.
1529
1531{
1532 SetBit(kMoveOpaque,set);
1533}
1534
1535////////////////////////////////////////////////////////////////////////////////
1536/// Paint canvas.
1537
1539{
1540 if (fCanvas)
1542}
1543
1544////////////////////////////////////////////////////////////////////////////////
1545/// Prepare for pick, call TPad::Pick() and when selected object
1546/// is different from previous then emit Picked() signal.
1547
1549{
1550 TObjLink *pickobj = nullptr;
1551
1552 fSelected = nullptr;
1553 fSelectedOpt = "";
1554 fSelectedPad = nullptr;
1555
1556 TPad *pad = Pick(px, py, pickobj);
1557 if (!pad) return nullptr;
1558
1559 if (!pickobj) {
1560 fSelected = pad;
1561 fSelectedOpt = "";
1562 } else {
1563 if (!fSelected) { // can be set via TCanvas::SetSelected()
1564 fSelected = pickobj->GetObject();
1565 fSelectedOpt = pickobj->GetOption();
1566 }
1567 }
1568 fSelectedPad = pad;
1569
1570 if (fSelected != prevSelObj)
1571 Picked(fSelectedPad, fSelected, fEvent); // emit signal
1572
1573 if ((fEvent == kButton1Down) || (fEvent == kButton2Down) || (fEvent == kButton3Down)) {
1577 Selected(fSelectedPad, fSelected, fEvent); // emit signal
1578 fSelectedX = px;
1579 fSelectedY = py;
1580 }
1581 }
1582 return pad;
1583}
1584
1585////////////////////////////////////////////////////////////////////////////////
1586/// Emit Picked() signal.
1587
1589{
1590 Longptr_t args[3];
1591
1592 args[0] = (Longptr_t) pad;
1593 args[1] = (Longptr_t) obj;
1594 args[2] = event;
1595
1596 Emit("Picked(TPad*,TObject*,Int_t)", args);
1597}
1598
1599////////////////////////////////////////////////////////////////////////////////
1600/// Emit Highlighted() signal.
1601///
1602/// - pad is pointer to pad with highlighted histogram or graph
1603/// - obj is pointer to highlighted histogram or graph
1604/// - x is highlighted x bin for 1D histogram or highlighted x-th point for graph
1605/// - y is highlighted y bin for 2D histogram (for 1D histogram or graph not in use)
1606
1608{
1609 Longptr_t args[4];
1610
1611 args[0] = (Longptr_t) pad;
1612 args[1] = (Longptr_t) obj;
1613 args[2] = x;
1614 args[3] = y;
1615
1616 Emit("Highlighted(TVirtualPad*,TObject*,Int_t,Int_t)", args);
1617}
1618
1619////////////////////////////////////////////////////////////////////////////////
1620/// This is "simplification" for function TCanvas::Connect with Highlighted
1621/// signal for specific slot.
1622///
1623/// Slot has to be defined "UserFunction(TVirtualPad *pad, TObject *obj, Int_t x, Int_t y)"
1624/// all parameters of UserFunction are taken from TCanvas::Highlighted
1625
1627{
1628 Connect("Highlighted(TVirtualPad*,TObject*,Int_t,Int_t)", nullptr, nullptr, slot);
1629}
1630
1631////////////////////////////////////////////////////////////////////////////////
1632/// Emit Selected() signal.
1633
1635{
1636 Longptr_t args[3];
1637
1638 args[0] = (Longptr_t) pad;
1639 args[1] = (Longptr_t) obj;
1640 args[2] = event;
1641
1642 Emit("Selected(TVirtualPad*,TObject*,Int_t)", args);
1643}
1644
1645////////////////////////////////////////////////////////////////////////////////
1646/// Emit ProcessedEvent() signal.
1647
1649{
1650 Longptr_t args[4];
1651
1652 args[0] = event;
1653 args[1] = x;
1654 args[2] = y;
1655 args[3] = (Longptr_t) obj;
1656
1657 Emit("ProcessedEvent(Int_t,Int_t,Int_t,TObject*)", args);
1658}
1659
1660////////////////////////////////////////////////////////////////////////////////
1661/// Recompute canvas parameters following a X11 Resize.
1662
1664{
1665 if (fCanvasID == -1)
1666 return;
1667
1668 if (!gROOT->IsLineProcessing() && !gVirtualX->IsCmdThread()) {
1669 gInterpreter->Execute(this, IsA(), "Resize", "");
1670 return;
1671 }
1672
1674
1675 TContext ctxt(this, kTRUE);
1676
1677 if (!IsBatch() && !IsWeb()) {
1678 // SL: do we need it here?
1679 fPainter->SelectDrawable(fCanvasID); //select current canvas for painting???
1680
1681 fCanvasImp->ResizeCanvasWindow(fCanvasID); //resize canvas and off-screen buffer
1682
1683 // Get effective window parameters including menubar and borders
1686
1687 // Get effective canvas parameters without borders
1689 }
1690
1691 if (fXsizeUser && fYsizeUser) {
1692 UInt_t nwh = fCh;
1693 UInt_t nww = fCw;
1695 if (rxy < 1) {
1697 if (twh > fCh)
1699 else
1700 nwh = twh;
1701 if (nww > fCw) {
1702 nww = fCw; nwh = twh;
1703 }
1704 if (nwh > fCh) {
1705 nwh = fCh; nww = UInt_t(Double_t(fCh)/rxy);
1706 }
1707 } else {
1709 if (twh > fCh)
1711 else
1712 nww = twh;
1713 if (nww > fCw) {
1714 nww = fCw; nwh = twh;
1715 }
1716 if (nwh > fCh) {
1717 nwh = fCh; nww = UInt_t(Double_t(fCh)*rxy);
1718 }
1719 }
1720 fCw = nww;
1721 fCh = nwh;
1722 }
1723
1724 if (fCw < fCh) {
1727 }
1728 else {
1731 }
1732
1733 //*-*- Loop on all pads to recompute conversion coefficients
1735}
1736
1737
1738////////////////////////////////////////////////////////////////////////////////
1739/// Raise canvas window
1740
1742{
1743 if (fCanvasImp)
1745}
1746
1747////////////////////////////////////////////////////////////////////////////////
1748/// Set option to resize objects/pads in a canvas.
1749///
1750/// - set = 1 (default) graphics objects are resized in opaque mode
1751/// - set = 0 only the outline of objects is drawn when resizing them
1752///
1753/// The option opaque produces the best effect. It requires however a
1754/// a reasonably fast workstation or response time.
1755
1757{
1758 SetBit(kResizeOpaque,set);
1759}
1760
1761////////////////////////////////////////////////////////////////////////////////
1762/// Execute the list of TExecs in the current pad.
1763
1765{
1766 if (!TestBit(kAutoExec))
1767 return;
1768 if (gPad)
1769 ((TPad*)gPad)->AutoExec();
1770}
1771
1772////////////////////////////////////////////////////////////////////////////////
1773/// Save primitives in this canvas in C++ macro file with GUI.
1774
1775void TCanvas::SavePrimitive(std::ostream &out, Option_t *option /*= ""*/)
1776{
1777 // Write canvas options (in $TROOT or $TStyle)
1778 out << " gStyle->SetOptFit(" << gStyle->GetOptFit() << ");\n";
1779 out << " gStyle->SetOptStat(" << gStyle->GetOptStat() << ");\n";
1780 out << " gStyle->SetOptTitle(" << gStyle->GetOptTitle() << ");\n";
1781
1782 if (gROOT->GetEditHistograms())
1783 out << " gROOT->SetEditHistograms();\n";
1784
1785 if (GetShowEventStatus())
1786 out << " " << GetName() << "->ToggleEventStatus();\n";
1787
1788 if (GetShowToolTips())
1789 out << " " << GetName() << "->ToggleToolTips();\n";
1790
1791 if (GetShowToolBar())
1792 out << " " << GetName() << "->ToggleToolBar();\n";
1793 if (GetHighLightColor() != 5)
1794 out << " " << GetName() << "->SetHighLightColor(" << TColor::SavePrimitiveColor(GetHighLightColor()) << ");\n";
1795
1796 // Now recursively scan all pads of this canvas
1797 cd();
1799}
1800
1801////////////////////////////////////////////////////////////////////////////////
1802/// Save primitives in this canvas as a C++ macro file.
1803/// This function loops on all the canvas primitives and for each primitive
1804/// calls the object SavePrimitive function.
1805/// When outputting floating point numbers, the default precision is 7 digits.
1806/// The precision can be changed (via system.rootrc) by changing the value
1807/// of the environment variable "Canvas.SavePrecision"
1808
1809void TCanvas::SaveSource(const char *filename, Option_t * /*option*/)
1810{
1811 // Reset the ClassSaved status of all classes
1812 gROOT->ResetClassSaved();
1813
1816
1818 if (cname.IsNull()) {
1819 invalid = kTRUE;
1820 cname = "c1";
1821 }
1822
1823 // if filename is given, open this file, otherwise create a file
1824 // with a name equal to the canvasname.C
1825 TString fname;
1826 if (filename && *filename) {
1827 fname = filename;
1828 } else {
1829 fname = cname + ".C";
1830 }
1831
1832 std::ofstream out;
1833 out.open(fname.Data(), std::ios::out);
1834 if (!out.good()) {
1835 Error("SaveSource", "Cannot open file: %s", fname.Data());
1836 return;
1837 }
1838
1839 //set precision
1840 Int_t precision = gEnv->GetValue("Canvas.SavePrecision",7);
1841 out.precision(precision);
1842
1843 // Write macro header and date/time stamp
1844 TDatime t;
1846 Int_t topx,topy;
1847 UInt_t w, h;
1848 if (!fCanvasImp) {
1849 Error("SaveSource", "Cannot open TCanvas");
1850 return;
1851 }
1854 h = UInt_t((fWindowHeight)/cx);
1855 topx = GetWindowTopX();
1856 topy = GetWindowTopY();
1857
1858 if (w == 0) {
1859 w = GetWw()+4; h = GetWh()+4;
1860 topx = 1; topy = 1;
1861 }
1862
1864 out << R"CODE(#ifdef __CLING__
1865#pragma cling optimize(0)
1866#endif
1867)CODE";
1868 Int_t p = mname.Last('.');
1869 Int_t s = mname.Last('/')+1;
1870
1871 // A named macro is generated only if the function name is valid. If not, the
1872 // macro is unnamed.
1873 TString first(mname(s,s+1));
1874 if (!first.IsDigit()) out <<"void " << mname(s,p-s) << "()" << std::endl;
1875
1876 out <<"{"<<std::endl;
1877 out <<"//=========Macro generated from canvas: "<<GetName()<<"/"<<GetTitle()<<std::endl;
1878 out <<"//========= ("<<t.AsString()<<") by ROOT version "<<gROOT->GetVersion()<<std::endl;
1879
1881 out <<std::endl<<" gStyle->SetCanvasPreferGL(kTRUE);"<<std::endl<<std::endl;
1882
1883 // Write canvas parameters (TDialogCanvas case)
1885 out << " " << ClassName() << " *" << cname << " = new " << ClassName() << "(\"" << GetName() << "\", \""
1886 << TString(GetTitle()).ReplaceSpecialCppChars() << "\", " << w << ", " << h << ");\n";
1887 } else {
1888 // Write canvas parameters (TCanvas case)
1889 out << " TCanvas *" << cname << " = new TCanvas(\"" << GetName() << "\", \""
1890 << TString(GetTitle()).ReplaceSpecialCppChars() << "\", " << (HasMenuBar() ? topx : -topx) << ", " << topy
1891 << ", " << w << ", " << h << ");\n";
1892 }
1893 // Write canvas options (in $TROOT or $TStyle)
1894 out << " gStyle->SetOptFit(" << gStyle->GetOptFit() << ");\n";
1895 out << " gStyle->SetOptStat(" << gStyle->GetOptStat() << ");\n";
1896 out << " gStyle->SetOptTitle(" << gStyle->GetOptTitle() << ");\n";
1897 if (gROOT->GetEditHistograms())
1898 out << " gROOT->SetEditHistograms();\n";
1899 if (GetShowEventStatus())
1900 out << " " << GetName() << "->ToggleEventStatus();\n";
1901 if (GetShowToolTips())
1902 out << " " << GetName() << "->ToggleToolTips();\n";
1903 if (GetHighLightColor() != 5)
1904 out << " " << GetName() << "->SetHighLightColor(" << TColor::SavePrimitiveColor(GetHighLightColor()) << ");\n";
1905
1907
1908 // Now recursively scan all pads of this canvas
1909 cd();
1910 if (invalid) fName = cname;
1911 TPad::SavePrimitive(out,"toplevel");
1912
1913 // Write canvas options related to pad editor
1914 out << " " << GetName() << "->SetSelected(" << GetName() << ");\n";
1915 if (GetShowToolBar())
1916 out << " " << GetName() << "->ToggleToolBar();\n";
1917 if (invalid)
1918 fName = cname0;
1919
1920 out <<"}\n";
1921 out.close();
1922 Info("SaveSource","C++ Macro file: %s has been generated", fname.Data());
1923
1924 // Reset the ClassSaved status of all classes
1925 gROOT->ResetClassSaved();
1926}
1927
1928////////////////////////////////////////////////////////////////////////////////
1929/// Toggle batch mode. However, if the canvas is created without a window
1930/// then batch mode always stays set.
1931
1933{
1934 if (gROOT->IsBatch() || IsWeb())
1935 fBatch = kTRUE;
1936 else
1937 fBatch = batch;
1938}
1939
1940////////////////////////////////////////////////////////////////////////////////
1941/// Set Width and Height of canvas to ww and wh respectively. If ww and/or wh
1942/// are greater than the current canvas window a scroll bar is automatically
1943/// generated. Use this function to zoom in a canvas and navigate via
1944/// the scroll bars. The Width and Height in this method are different from those
1945/// given in the TCanvas constructors where these two dimension include the size
1946/// of the window decoration whereas they do not in this method.
1947/// When both ww==0 and wh==0, auto resize mode will be enabled again and
1948/// canvas drawing area will automatically fit available window size
1949
1951{
1952 if (fCanvasImp) {
1953 fCw = ww;
1954 fCh = wh;
1956 TContext ctxt(this, kTRUE);
1957 ResizePad();
1958 }
1959}
1960
1961////////////////////////////////////////////////////////////////////////////////
1962/// Set cursor.
1963
1965{
1966 if (!IsBatch() && !IsWeb() && fCanvasID != -1)
1968}
1969
1970////////////////////////////////////////////////////////////////////////////////
1971/// Set Double Buffer On/Off.
1972
1974{
1975 if (IsBatch() || IsWeb())
1976 return;
1978 if (fCanvasID != -1)
1980
1981 // depending of the buffer mode set the drawing window to either
1982 // the canvas pixmap or to the canvas on-screen window
1983 if (fDoubleBuffer) {
1985 } else
1987}
1988
1989////////////////////////////////////////////////////////////////////////////////
1990/// Fix canvas aspect ratio to current value if fixed is true.
1991
1993{
1994 if (fixed) {
1995 if (!fFixedAspectRatio) {
1996 if (fCh != 0)
1998 else {
1999 Error("SetAspectRatio", "cannot fix aspect ratio, height of canvas is 0");
2000 return;
2001 }
2003 }
2004 } else {
2006 fAspectRatio = 0;
2007 }
2008}
2009
2010////////////////////////////////////////////////////////////////////////////////
2011/// If isfolder=kTRUE, the canvas can be browsed like a folder
2012/// by default a canvas is not browsable.
2013
2015{
2017}
2018
2019////////////////////////////////////////////////////////////////////////////////
2020/// Set canvas name. In case `name` is an empty string, a default name is set.
2021/// Canvas automatically marked as modified when SetName method called
2022
2023void TCanvas::SetName(const char *name)
2024{
2026
2027 Modified();
2028}
2029
2030
2031////////////////////////////////////////////////////////////////////////////////
2032/// Function to resize a canvas so that the plot inside is shown in real aspect
2033/// ratio
2034///
2035/// \param[in] axis 1 for resizing horizontally (x-axis) in order to get real
2036/// aspect ratio, 2 for the resizing vertically (y-axis)
2037/// \return false if error is encountered, true otherwise
2038///
2039/// ~~~ {.cpp}
2040/// hpxpy->Draw();
2041/// c1->SetRealAspectRatio();
2042/// ~~~
2043///
2044/// - For defining the concept of real aspect ratio, it is assumed that x and y
2045/// axes are in same units, e.g. both in MeV or both in ns.
2046/// - You can resize either the width of the canvas or the height, but not both
2047/// at the same time
2048/// - Call this function AFTER drawing AND zooming (SetUserRange) your TGraph or
2049/// Histogram, otherwise it cannot infer your actual axes lengths
2050/// - This function ensures that the TFrame has a real aspect ratio, this does not
2051/// mean that the full pad (i.e. the canvas or png output) including margins has
2052/// exactly the same ratio
2053/// - This function does not work if the canvas is divided in several subpads
2054
2055bool TCanvas::SetRealAspectRatio(const Int_t axis)
2056{
2057 Update();
2058
2059 //Get how many pixels are occupied by the canvas
2060 Int_t npx = GetWw();
2061 Int_t npy = GetWh();
2062
2063 //Get x-y coordinates at the edges of the canvas (extrapolating outside the axes, NOT at the edges of the histogram)
2064 Double_t x1 = GetX1();
2065 Double_t y1 = GetY1();
2066 Double_t x2 = GetX2();
2067 Double_t y2 = GetY2();
2068
2069 //Get the length of extrapolated x and y axes
2070 Double_t xlength2 = x2 - x1;
2071 Double_t ylength2 = y2 - y1;
2073
2074 //Now get the number of pixels including the canvas borders
2077
2078 if (axis==1) {
2081 } else if (axis==2) {
2084 } else {
2085 Error("SetRealAspectRatio", "axis value %d is neither 1 (resize along x-axis) nor 2 (resize along y-axis).",axis);
2086 return false;
2087 }
2088
2089 //Check now that resizing has worked
2090
2091 Update();
2092
2093 //Get how many pixels are occupied by the canvas
2094 npx = GetWw();
2095 npy = GetWh();
2096
2097 //Get x-y coordinates at the edges of the canvas (extrapolating outside the axes,
2098 //NOT at the edges of the histogram)
2099 x1 = GetX1();
2100 y1 = GetY1();
2101 x2 = GetX2();
2102 y2 = GetY2();
2103
2104 //Get the length of extrapolated x and y axes
2105 xlength2 = x2 - x1;
2106 ylength2 = y2 - y1;
2108
2109 //Check accuracy +/-1 pixel due to rounding
2110 if (abs(TMath::Nint(npy*ratio2) - npx)<2) {
2111 return true;
2112 } else {
2113 Error("SetRealAspectRatio", "Resizing failed.");
2114 return false;
2115 }
2116}
2117
2118
2119////////////////////////////////////////////////////////////////////////////////
2120/// Set selected canvas.
2121
2123{
2124 fSelected = obj;
2125 if (obj) obj->SetBit(kMustCleanup);
2126}
2127
2128////////////////////////////////////////////////////////////////////////////////
2129/// Set canvas title.
2130
2131void TCanvas::SetTitle(const char *title)
2132{
2133 fTitle = title;
2135}
2136
2137////////////////////////////////////////////////////////////////////////////////
2138/// Set canvas window position
2139
2141{
2142 if (fCanvasImp)
2144}
2145
2146////////////////////////////////////////////////////////////////////////////////
2147/// Set canvas window size
2148
2150{
2151 if (fBatch && !IsWeb())
2152 SetCanvasSize((ww + fCw) / 2, (wh + fCh) / 2);
2153 else if (fCanvasImp)
2155}
2156
2157////////////////////////////////////////////////////////////////////////////////
2158/// Set canvas implementation
2159/// If web-based implementation provided, some internal fields also initialized
2160
2162{
2163 Bool_t was_web = IsWeb();
2164
2165 fCanvasImp = imp;
2166
2167 if (!was_web && IsWeb()) {
2169 fPixmapID = 0;
2170 fMother = this;
2171 if (!fCw) fCw = 800;
2172 if (!fCh) fCh = 600;
2173 } else if (was_web && !imp) {
2174 fCanvasID = -1;
2175 fPixmapID = -1;
2176 fMother = nullptr;
2177 fCw = fCh = 0;
2178 }
2179}
2180
2181////////////////////////////////////////////////////////////////////////////////
2182/// Set the canvas scale in centimeters.
2183///
2184/// This information is used by PostScript to set the page size.
2185///
2186/// \param[in] xsize size of the canvas in centimeters along X
2187/// \param[in] ysize size of the canvas in centimeters along Y
2188///
2189/// if xsize and ysize are not equal to 0, then the scale factors will
2190/// be computed to keep the ratio ysize/xsize independently of the canvas
2191/// size (parts of the physical canvas will be unused).
2192///
2193/// if xsize = 0 and ysize is not zero, then xsize will be computed
2194/// to fit to the current canvas scale. If the canvas is resized,
2195/// a new value for xsize will be recomputed. In this case the aspect
2196/// ratio is not preserved.
2197///
2198/// if both xsize = 0 and ysize = 0, then the scaling is automatic.
2199/// the largest dimension will be allocated a size of 20 centimeters.
2200
2202{
2203 fXsizeUser = xsize;
2204 fYsizeUser = ysize;
2205
2206 Resize();
2207}
2208
2209////////////////////////////////////////////////////////////////////////////////
2210/// Show canvas
2211
2212void TCanvas::Show()
2213{
2214 if (fCanvasImp)
2215 fCanvasImp->Show();
2216}
2217
2218////////////////////////////////////////////////////////////////////////////////
2219/// Stream a class object.
2220
2222{
2223 UInt_t R__s, R__c;
2224 if (b.IsReading()) {
2225 Version_t v = b.ReadVersion(&R__s, &R__c);
2226 gPad = this;
2227 fCanvas = this;
2228 if (v>7) b.ClassBegin(TCanvas::IsA());
2229 if (v>7) b.ClassMember("TPad");
2231 gPad = this;
2232 //restore the colors
2233 auto colors = dynamic_cast<TObjArray *>(fPrimitives->FindObject("ListOfColors"));
2234 if (colors) {
2235 auto root_colors = dynamic_cast<TObjArray *>(gROOT->GetListOfColors());
2236
2237 TIter next(colors);
2238 while (auto colold = static_cast<TColor *>(next())) {
2239 Int_t cn = colold->GetNumber();
2240 TColor *colcur = gROOT->GetColor(cn);
2241 if (colcur && (colcur->IsA() == TColor::Class()) && (colold->IsA() == TColor::Class())) {
2242 colcur->SetName(colold->GetName());
2243 colcur->SetRGB(colold->GetRed(), colold->GetGreen(), colold->GetBlue());
2244 colcur->SetAlpha(colold->GetAlpha());
2245 } else {
2246 if (colcur) {
2247 if (root_colors) root_colors->Remove(colcur);
2248 delete colcur;
2249 }
2250 colors->Remove(colold);
2251 if (root_colors) {
2252 if (colcur) {
2253 root_colors->AddAtAndExpand(colold, cn);
2254 }
2255 else {
2256 // Copy to current session
2257 // do not use copy constructor which does not update highest color index
2258 [[maybe_unused]] TColor* const colnew = new TColor(cn, colold->GetRed(), colold->GetGreen(), colold->GetBlue(), colold->GetName(), colold->GetAlpha());
2259 delete colold;
2260 // No need to delete colnew, as the constructor adds it to global list of colors
2261 assert(root_colors->At(cn) == colnew);
2262 }
2263 }
2264 }
2265 }
2266 //restore the palette if needed
2267 auto palette = dynamic_cast<TObjArray *>(fPrimitives->FindObject("CurrentColorPalette"));
2268 if (palette) {
2270 Int_t number = palette->GetEntries();
2271 TArrayI palcolors(number);
2272 Int_t i = 0;
2273 while (auto col = static_cast<TColor *>(nextcol()))
2274 palcolors[i++] = col->GetNumber();
2275 gStyle->SetPalette(number, palcolors.GetArray());
2277 delete palette;
2278 }
2280 colors->Delete();
2281 delete colors;
2282 }
2283
2284 if (v>7) b.ClassMember("fDISPLAY","TString");
2286 if (v>7) b.ClassMember("fDoubleBuffer", "Int_t");
2287 b >> fDoubleBuffer;
2288 if (v>7) b.ClassMember("fRetained", "Bool_t");
2289 b >> fRetained;
2290 if (v>7) b.ClassMember("fXsizeUser", "Size_t");
2291 b >> fXsizeUser;
2292 if (v>7) b.ClassMember("fYsizeUser", "Size_t");
2293 b >> fYsizeUser;
2294 if (v>7) b.ClassMember("fXsizeReal", "Size_t");
2295 b >> fXsizeReal;
2296 if (v>7) b.ClassMember("fYsizeReal", "Size_t");
2297 b >> fYsizeReal;
2298 fCanvasID = -1;
2299 if (v>7) b.ClassMember("fWindowTopX", "Int_t");
2300 b >> fWindowTopX;
2301 if (v>7) b.ClassMember("fWindowTopY", "Int_t");
2302 b >> fWindowTopY;
2303 if (v > 2) {
2304 if (v>7) b.ClassMember("fWindowWidth", "UInt_t");
2305 b >> fWindowWidth;
2306 if (v>7) b.ClassMember("fWindowHeight", "UInt_t");
2307 b >> fWindowHeight;
2308 }
2309 if (v>7) b.ClassMember("fCw", "UInt_t");
2310 b >> fCw;
2311 if (v>7) b.ClassMember("fCh", "UInt_t");
2312 b >> fCh;
2313 if (v <= 2) {
2314 fWindowWidth = fCw;
2316 }
2317 if (v>7) b.ClassMember("fCatt", "TAttCanvas");
2318 fCatt.Streamer(b);
2319 Bool_t dummy;
2320 if (v>7) b.ClassMember("kMoveOpaque", "Bool_t");
2321 b >> dummy; if (dummy) MoveOpaque(1);
2322 if (v>7) b.ClassMember("kResizeOpaque", "Bool_t");
2323 b >> dummy; if (dummy) ResizeOpaque(1);
2324 if (v>7) b.ClassMember("fHighLightColor", "Color_t");
2325 b >> fHighLightColor;
2326 if (v>7) b.ClassMember("fBatch", "Bool_t");
2327 b >> dummy; //was fBatch
2328 if (v < 2) return;
2329 if (v>7) b.ClassMember("kShowEventStatus", "Bool_t");
2330 b >> dummy; if (dummy) SetBit(kShowEventStatus);
2331
2332 if (v > 3) {
2333 if (v>7) b.ClassMember("kAutoExec", "Bool_t");
2334 b >> dummy; if (dummy) SetBit(kAutoExec);
2335 }
2336 if (v>7) b.ClassMember("kMenuBar", "Bool_t");
2337 b >> dummy; if (dummy) SetBit(kMenuBar);
2338 fBatch = gROOT->IsBatch();
2339 if (v>7) b.ClassEnd(TCanvas::IsA());
2340 b.CheckByteCount(R__s, R__c, TCanvas::IsA());
2341 } else {
2342 //save list of colors
2343 //we must protect the case when two or more canvases are saved
2344 //in the same buffer. If the list of colors has already been saved
2345 //in the buffer, do not add the list of colors to the list of primitives.
2346 TObjArray *colors = nullptr;
2347 TObjArray *CurrentColorPalette = nullptr;
2348 if (TColor::DefinedColors()) {
2349 if (!b.CheckObject(gROOT->GetListOfColors(),TObjArray::Class())) {
2350 colors = (TObjArray*)gROOT->GetListOfColors();
2352 }
2353 //save the current palette
2355 Int_t palsize = pal.GetSize();
2357 CurrentColorPalette->SetName("CurrentColorPalette");
2358 for (Int_t i=0; i<palsize; i++) CurrentColorPalette->Add(gROOT->GetColor(pal[i]));
2360 }
2361
2362 R__c = b.WriteVersion(TCanvas::IsA(), kTRUE);
2363 b.ClassBegin(TCanvas::IsA());
2364 b.ClassMember("TPad");
2368 b.ClassMember("fDISPLAY","TString");
2370 b.ClassMember("fDoubleBuffer", "Int_t");
2371 b << fDoubleBuffer;
2372 b.ClassMember("fRetained", "Bool_t");
2373 b << fRetained;
2374 b.ClassMember("fXsizeUser", "Size_t");
2375 b << fXsizeUser;
2376 b.ClassMember("fYsizeUser", "Size_t");
2377 b << fYsizeUser;
2378 b.ClassMember("fXsizeReal", "Size_t");
2379 b << fXsizeReal;
2380 b.ClassMember("fYsizeReal", "Size_t");
2381 b << fYsizeReal;
2384 UInt_t editorWidth = 0;
2386 b.ClassMember("fWindowTopX", "Int_t");
2387 b << topx;
2388 b.ClassMember("fWindowTopY", "Int_t");
2389 b << topy;
2390 b.ClassMember("fWindowWidth", "UInt_t");
2391 b << (UInt_t)(w-editorWidth);
2392 b.ClassMember("fWindowHeight", "UInt_t");
2393 b << h;
2394 b.ClassMember("fCw", "UInt_t");
2395 b << fCw;
2396 b.ClassMember("fCh", "UInt_t");
2397 b << fCh;
2398 b.ClassMember("fCatt", "TAttCanvas");
2399 fCatt.Streamer(b);
2400 b.ClassMember("kMoveOpaque", "Bool_t");
2401 b << TestBit(kMoveOpaque); //please remove in ROOT version 6
2402 b.ClassMember("kResizeOpaque", "Bool_t");
2403 b << TestBit(kResizeOpaque); //please remove in ROOT version 6
2404 b.ClassMember("fHighLightColor", "Color_t");
2405 b << fHighLightColor;
2406 b.ClassMember("fBatch", "Bool_t");
2407 b << fBatch; //please remove in ROOT version 6
2408 b.ClassMember("kShowEventStatus", "Bool_t");
2409 b << TestBit(kShowEventStatus); //please remove in ROOT version 6
2410 b.ClassMember("kAutoExec", "Bool_t");
2411 b << TestBit(kAutoExec); //please remove in ROOT version 6
2412 b.ClassMember("kMenuBar", "Bool_t");
2413 b << TestBit(kMenuBar); //please remove in ROOT version 6
2414 b.ClassEnd(TCanvas::IsA());
2415 b.SetByteCount(R__c, kTRUE);
2416 }
2417}
2418
2419////////////////////////////////////////////////////////////////////////////////
2420/// Toggle pad auto execution of list of TExecs.
2421
2423{
2426}
2427
2428////////////////////////////////////////////////////////////////////////////////
2429/// Toggle event statusbar.
2430
2432{
2435
2437}
2438
2439////////////////////////////////////////////////////////////////////////////////
2440/// Toggle toolbar.
2441
2443{
2446
2448}
2449
2450////////////////////////////////////////////////////////////////////////////////
2451/// Toggle editor.
2452
2454{
2457
2459}
2460
2461////////////////////////////////////////////////////////////////////////////////
2462/// Toggle tooltip display.
2463
2465{
2468
2470}
2471
2472
2473////////////////////////////////////////////////////////////////////////////////
2474/// Static function returning "true" if transparency is supported.
2475
2477{
2478 if (gPad)
2479 if (auto pp = gPad->GetPainter())
2480 return pp->IsSupportAlpha();
2481 return kFALSE;
2482}
2483
2484extern "C" void ROOT_TCanvas_Update(void* TheCanvas) {
2485 static_cast<TCanvas*>(TheCanvas)->Update();
2486}
2487
2488////////////////////////////////////////////////////////////////////////////////
2489/// Update canvas pad buffers.
2490
2491void TCanvas::Update()
2492{
2493 fUpdated = kTRUE;
2494
2495 if (fUpdating) return;
2496
2497 if (fPixmapID == -1) return;
2498
2499 static const union CastFromFuncToVoidPtr_t {
2501 void (*fFuncPtr)(void*);
2502 void* fVoidPtr;
2504
2505 if (gThreadXAR) {
2506 void *arr[3];
2507 arr[1] = this;
2508 arr[2] = castFromFuncToVoidPtr.fVoidPtr;
2509 if ((*gThreadXAR)("CUPD", 3, arr, nullptr)) return;
2510 }
2511
2512 if (!fCanvasImp) return;
2513
2514 if (!gVirtualX->IsCmdThread()) {
2515 // Why do we have this (which uses the interpreter to funnel the Update()
2516 // through the main thread) when the gThreadXAR mechanism does seemingly
2517 // the same?
2518 gInterpreter->Execute(this, IsA(), "Update", "");
2519 return;
2520 }
2521
2523
2524 fUpdating = kTRUE;
2525
2527
2528 if (!IsBatch()) FeedbackMode(kFALSE); // Goto double buffer mode
2529
2530 if (!UseGL() || fGLDevice == -1) PaintModified(); // Repaint all modified pad's
2531
2532 Flush(); // Copy all pad pixmaps to the screen
2533
2535 }
2536
2537 fUpdating = kFALSE;
2538}
2539
2540////////////////////////////////////////////////////////////////////////////////
2541/// Asynchronous pad update.
2542/// In case of web-based canvas triggers update of the canvas on the client side,
2543/// but does not wait that real update is completed. Avoids blocking of caller thread.
2544/// Have to be used if called from other web-based widget to avoid logical dead-locks.
2545/// In case of normal canvas just canvas->Update() is performed.
2546
2548{
2549 fUpdated = kTRUE;
2550
2551 if (IsWeb())
2553 else
2554 Update();
2555}
2556
2557////////////////////////////////////////////////////////////////////////////////
2558/// Used by friend class TCanvasImp.
2559
2561{
2562 fCanvasID = 0;
2563 fContextMenu = nullptr;
2564}
2565
2566////////////////////////////////////////////////////////////////////////////////
2567/// Check whether this canvas is to be drawn in grayscale mode.
2568
2570{
2571 return TestBit(kIsGrayscale);
2572}
2573
2574////////////////////////////////////////////////////////////////////////////////
2575/// Set whether this canvas should be painted in grayscale, and re-paint
2576/// it if necessary.
2577
2578void TCanvas::SetGrayscale(Bool_t set /*= kTRUE*/)
2579{
2580 if (IsGrayscale() == set) return;
2581 SetBit(kIsGrayscale, set);
2582 if (IsWeb()) {
2583 Modified();
2584 UpdateAsync();
2585 } else {
2586 Paint(); // update canvas and all sub-pads, unconditionally!
2587 }
2588}
2589
2590////////////////////////////////////////////////////////////////////////////////
2591/// Probably, TPadPainter must be placed in a separate ROOT module -
2592/// "padpainter" (the same as "histpainter"). But now, it's directly in a
2593/// gpad dir, so, in case of default painter, no *.so should be loaded,
2594/// no need in plugin managers.
2595/// May change in future.
2596
2598{
2599 //Even for batch mode painter is still required, just to delegate
2600 //some calls to batch "virtual X".
2601 if (!UseGL() || fBatch || IsWeb()) {
2602 fPainter = nullptr;
2604 if (!fPainter) fPainter = new TPadPainter; // Do not need plugin manager for this!
2605 } else {
2607 if (!fPainter) {
2608 Error("CreatePainter", "GL Painter creation failed! Will use default!");
2609 fPainter = new TPadPainter;
2610 fUseGL = kFALSE;
2611 }
2612 }
2613}
2614
2615////////////////////////////////////////////////////////////////////////////////
2616/// Access and (probably) creation of pad painter.
2617
2619{
2620 if (!fPainter) CreatePainter();
2621 return fPainter;
2622}
2623
2624////////////////////////////////////////////////////////////////////////////////
2625/// Replace canvas painter
2626/// For intenral use only - when creating PS images
2627
2629{
2630 if (!create) {
2631 delete fPainter;
2632 fPainter = oldp;
2633 return kFALSE;
2634 }
2635
2636 if (!gVirtualPS /* || !IsBatch() */)
2637 return kFALSE;
2638
2639
2641 return kFALSE;
2642
2643 oldp = fPainter;
2645 return kTRUE;
2646}
2647
2648////////////////////////////////////////////////////////////////////////////////
2649///assert on IsBatch() == false?
2650
2652{
2653 if (fGLDevice != -1) {
2654 //fPainter has a font manager.
2655 //Font manager will delete textures.
2656 //If context is wrong (we can have several canvases) -
2657 //wrong texture will be deleted, damaging some of our fonts.
2658 gGLManager->MakeCurrent(fGLDevice);
2659 }
2660
2662
2663 if (fGLDevice != -1) {
2664 gGLManager->DeleteGLContext(fGLDevice);//?
2665 fGLDevice = -1;
2666 }
2667}
2668
2669
2670////////////////////////////////////////////////////////////////////////////////
2671/// Save provided pads/canvases into the image file(s)
2672/// Filename can include printf argument for image number - like "image%03d.png".
2673/// In this case images: "image000.png", "image001.png", "image002.png" will be created.
2674/// If pattern is not provided - it will be automatically inserted before extension except PDF and ROOT files.
2675/// In last case PDF or ROOT file will contain all pads.
2676/// Parameter option only used when output into PDF/PS files
2677/// If TCanvas::SaveAll() called without arguments - all existing canvases will be stored in allcanvases.pdf file.
2678
2679Bool_t TCanvas::SaveAll(const std::vector<TPad *> &pads, const char *filename, Option_t *option)
2680{
2681 if (pads.empty()) {
2682 std::vector<TPad *> canvases;
2683 TIter iter(gROOT->GetListOfCanvases());
2684 while (auto c = dynamic_cast<TCanvas *>(iter()))
2685 canvases.emplace_back(c);
2686
2687 if (canvases.empty()) {
2688 ::Warning("TCanvas::SaveAll", "No pads are provided");
2689 return kFALSE;
2690 }
2691
2692 return TCanvas::SaveAll(canvases, filename && *filename ? filename : "allcanvases.pdf", option);
2693 }
2694
2696
2697 Bool_t hasArg = fname.Contains("%");
2698
2699 if ((pads.size() == 1) && !hasArg) {
2700 pads[0]->SaveAs(filename);
2701 return kTRUE;
2702 }
2703
2704 auto p = fname.Last('.');
2705 if (p != kNPOS) {
2706 ext = fname(p+1, fname.Length() - p - 1);
2707 ext.ToLower();
2708 } else {
2709 p = fname.Length();
2710 ::Warning("TCanvas::SaveAll", "Extension is not provided in file name %s, append .png", filename);
2711 fname.Append(".png");
2712 ext = "png";
2713 }
2714
2715 if (ext != "pdf" && ext != "ps" && ext != "root" && ext != "xml" && !hasArg) {
2716 fname.Insert(p, "%d");
2717 hasArg = kTRUE;
2718 }
2719
2720 static std::vector<TString> webExtensions = { "png", "json", "svg", "pdf", "jpg", "jpeg", "webp" };
2721
2722 if (gROOT->IsWebDisplay()) {
2724 for (auto &wext : webExtensions) {
2725 if ((isSupported = (wext == ext)))
2726 break;
2727 }
2728
2729 if (isSupported) {
2730 auto cmd = TString::Format("TWebCanvas::ProduceImages( *((std::vector<TPad *> *) 0x%zx), \"%s\")", (size_t) &pads, fname.Data());
2731
2732 return (Bool_t) gROOT->ProcessLine(cmd);
2733 }
2734
2735 if ((ext != "root") && (ext != "xml"))
2736 ::Warning("TCanvas::SaveAll", "TWebCanvas does not support image format %s - using normal ROOT functionality", fname.Data());
2737 }
2738
2739 // store all pads into single PDF/PS files
2740 if (ext == "pdf" || ext == "ps") {
2741 for (unsigned n = 0; n < pads.size(); ++n) {
2742 TString fn = fname;
2743 if (hasArg)
2744 fn = TString::Format(fname.Data(), (int) n);
2745 else if (n == 0)
2746 fn.Append("(");
2747 else if (n == pads.size() - 1)
2748 fn.Append(")");
2749
2750 pads[n]->Print(fn.Data(), option && *option ? option : ext.Data());
2751 }
2752
2753 return kTRUE;
2754 }
2755
2756 // store all pads in single ROOT file
2757 if ((ext == "root" || ext == "xml") && !hasArg) {
2758 TString fn = fname;
2760 if (fn.IsNull()) {
2761 fn.Form("%s.%s", pads[0]->GetName(), ext.Data());
2762 ::Warning("TCanvas::SaveAll", "Filename %s cannot be used - use pad name %s as pattern", fname.Data(), fn.Data());
2763 }
2764
2766
2767 if (!gDirectory) {
2768 isError = kTRUE;
2769 } else {
2770 for (unsigned n = 0; n < pads.size(); ++n) {
2771 auto sz = gDirectory->SaveObjectAs(pads[n], fn.Data(), n==0 ? "q" : "qa");
2772 if (!sz) { isError = kTRUE; break; }
2773 }
2774 }
2775
2776 if (isError)
2777 ::Error("TCanvas::SaveAll", "Failure to store pads in %s", filename);
2778 else
2779 ::Info("TCanvas::SaveAll", "ROOT file %s has been created", filename);
2780
2781 return !isError;
2782 }
2783
2784 for (unsigned n = 0; n < pads.size(); ++n) {
2785 TString fn = TString::Format(fname.Data(), (int) n);
2787 if (fn.IsNull()) {
2788 fn.Form("%s%d.%s", pads[n]->GetName(), (int) n, ext.Data());
2789 ::Warning("TCanvas::SaveAll", "Filename %s cannot be used - use pad name %s as pattern", fname.Data(), fn.Data());
2790 }
2791
2792 pads[n]->SaveAs(fn.Data());
2793 }
2794
2795 return kTRUE;
2796
2797}
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:373
@ kCross
Definition GuiTypes.h:375
#define SafeDelete(p)
Definition RConfig.hxx:531
#define b(i)
Definition RSha256.hxx:100
#define c(i)
Definition RSha256.hxx:101
#define h(i)
Definition RSha256.hxx:106
bool Bool_t
Boolean (0=false, 1=true) (bool)
Definition RtypesCore.h:77
int Int_t
Signed integer 4 bytes (int)
Definition RtypesCore.h:59
float Size_t
Attribute size (float)
Definition RtypesCore.h:103
long Longptr_t
Integer large enough to hold a pointer (platform-dependent)
Definition RtypesCore.h:89
short Version_t
Class version identifier (short)
Definition RtypesCore.h:79
unsigned int UInt_t
Unsigned integer 4 bytes (unsigned int)
Definition RtypesCore.h:60
float Float_t
Float 4 bytes (float)
Definition RtypesCore.h:71
constexpr Bool_t kFALSE
Definition RtypesCore.h:108
double Double_t
Double 8 bytes.
Definition RtypesCore.h:73
constexpr Ssiz_t kNPOS
The equivalent of std::string::npos for the ROOT class TString.
Definition RtypesCore.h:131
constexpr Bool_t kTRUE
Definition RtypesCore.h:107
const char Option_t
Option string (const char)
Definition RtypesCore.h:80
@ kRed
Definition Rtypes.h:67
R__EXTERN TApplication * gApplication
void ROOT_TCanvas_Update(void *TheCanvas)
Definition TCanvas.cxx:2481
class TCanvasInit gCanvasInit
TString GetNewCanvasName(const char *arg=nullptr)
Definition TCanvas.cxx:66
const Size_t kDefaultCanvasSize
Definition TCanvas.cxx:63
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
#define gDirectory
Definition TDirectory.h:385
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:157
R__EXTERN TGuiFactory * gBatchGuiFactory
Definition TGuiFactory.h:67
R__EXTERN TGuiFactory * gGuiFactory
Definition TGuiFactory.h:66
#define gInterpreter
Int_t gDebug
Global variable setting the debug level. Set to 0 to disable, increase it in steps of 1 to increase t...
Definition TROOT.cxx:627
R__EXTERN TVirtualMutex * gROOTMutex
Definition TROOT.h:63
#define gROOT
Definition TROOT.h:414
R__EXTERN TStyle * gStyle
Definition TStyle.h:442
R__EXTERN TSystem * gSystem
Definition TSystem.h:582
#define gGLManager
Definition TVirtualGL.h:159
#define R__LOCKGUARD(mutex)
R__EXTERN TVirtualPS * gVirtualPS
Definition TVirtualPS.h:84
#define gPad
R__EXTERN Int_t(* gThreadXAR)(const char *xact, Int_t nb, void **ar, Int_t *iret)
#define gVirtualX
Definition TVirtualX.h:365
Color * colors
Definition X3DBuffer.c:21
#define snprintf
Definition civetweb.c:1579
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
virtual void Streamer(TBuffer &)
Fill Area Attributes class.
Definition TAttFill.h:21
virtual Color_t GetFillColor() const
Return the fill area color.
Definition TAttFill.h:32
void Copy(TAttFill &attfill) const
Copy this fill attributes to a new TAttFill.
Definition TAttFill.cxx:206
Line Attributes class.
Definition TAttLine.h:21
void Copy(TAttLine &attline) const
Copy this line attributes to a new TAttLine.
Definition TAttLine.cxx:176
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:98
virtual void SetLeftMargin(Float_t leftmargin)
Set Pad left margin in fraction of the pad width.
Definition TAttPad.cxx:108
virtual void Copy(TAttPad &attpad) const
copy function
Definition TAttPad.cxx:43
virtual void SetRightMargin(Float_t rightmargin)
Set Pad right margin in fraction of the pad width.
Definition TAttPad.cxx:118
virtual void SetTopMargin(Float_t topmargin)
Set Pad top margin in fraction of the pad height.
Definition TAttPad.cxx:128
Class to manage histogram axis.
Definition TAxis.h:32
static TClass * Class()
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
ABC describing GUI independent main window (with menubar, scrollbars and a drawing area).
Definition TCanvasImp.h:30
virtual void SetStatusText(const char *text=nullptr, Int_t partidx=0)
Definition TCanvasImp.h:70
virtual void SetWindowPosition(Int_t x, Int_t y)
Definition TCanvasImp.h:71
virtual void ShowMenuBar(Bool_t show=kTRUE)
Definition TCanvasImp.h:76
virtual void Show()
Definition TCanvasImp.h:75
virtual Bool_t PerformUpdate(Bool_t)
Definition TCanvasImp.h:49
virtual void ShowToolTips(Bool_t show=kTRUE)
Definition TCanvasImp.h:83
virtual void Iconify()
Definition TCanvasImp.h:68
virtual Int_t InitWindow()
Definition TCanvasImp.h:69
virtual void UpdateDisplay(Int_t mode=0, Bool_t sleep=kFALSE)
Update gVirtualX display, also optionally sleep to wait until operation finished.
virtual void Close()
Definition TCanvasImp.h:59
virtual void SetWindowTitle(const char *newTitle)
Definition TCanvasImp.h:73
virtual void GetCanvasGeometry(Int_t wid, UInt_t &w, UInt_t &h)
Gets the size and position of the canvas paint area.
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 void ResizeCanvasWindow(Int_t wid)
Resize canvas window, redirect to gVirtualX.
virtual TVirtualPadPainter * CreatePadPainter()
Definition TCanvasImp.h:50
virtual void ShowEditor(Bool_t show=kTRUE)
Definition TCanvasImp.h:81
virtual void RaiseWindow()
Definition TCanvasImp.h:78
virtual void SetWindowSize(UInt_t width, UInt_t height)
Definition TCanvasImp.h:72
virtual void SetCanvasSize(UInt_t w, UInt_t h)
Definition TCanvasImp.h:74
virtual void ForceUpdate()
Definition TCanvasImp.h:60
virtual void ShowStatusBar(Bool_t show=kTRUE)
Definition TCanvasImp.h:77
virtual void ShowToolBar(Bool_t show=kTRUE)
Definition TCanvasImp.h:82
The Canvas class.
Definition TCanvas.h:23
void Init()
Initialize the TCanvas members. Called by all constructors.
Definition TCanvas.cxx:536
UInt_t fCw
Width of the canvas along X (pixels)
Definition TCanvas.h:44
~TCanvas() override
Canvas destructor.
Definition TCanvas.cxx:674
void EmbedInto(Int_t winid, Int_t ww, Int_t wh)
Embedded a canvas into a TRootEmbeddedCanvas.
Definition TCanvas.cxx:1047
void SetWindowSize(UInt_t ww, UInt_t wh)
Set canvas window size.
Definition TCanvas.cxx:2146
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:2011
void Browse(TBrowser *b) override
Browse.
Definition TCanvas.cxx:682
UInt_t GetWindowHeight() const
Definition TCanvas.h:164
virtual void EditorBar()
Get editor bar.
Definition TCanvas.cxx:1038
static TCanvas * MakeDefCanvas()
Static function to build a default canvas.
Definition TCanvas.cxx:1511
void EnterLeave(TPad *prevSelPad, TObject *prevSelObj)
Generate kMouseEnter and kMouseLeave events depending on the previously selected object and the curre...
Definition TCanvas.cxx:1073
Size_t fYsizeReal
Current size of canvas along Y in CM.
Definition TCanvas.h:37
void Constructor()
Canvas default constructor.
Definition TCanvas.cxx:190
virtual void ToggleAutoExec()
Toggle pad auto execution of list of TExecs.
Definition TCanvas.cxx:2419
TCanvas(const TCanvas &canvas)=delete
Int_t fWindowTopX
Top X position of window (in pixels)
Definition TCanvas.h:40
void Draw(Option_t *option="") override
Draw a canvas.
Definition TCanvas.cxx:853
void SetDoubleBuffer(Int_t mode=1) override
Set Double Buffer On/Off.
Definition TCanvas.cxx:1970
virtual void ToggleToolTips()
Toggle tooltip display.
Definition TCanvas.cxx:2461
void Clear(Option_t *option="") override
Remove all primitives from the canvas.
Definition TCanvas.cxx:734
void UseCurrentStyle() override
Force a copy of current style for all objects in canvas.
Definition TCanvas.cxx:1181
void Iconify()
Iconify canvas.
Definition TCanvas.cxx:1474
Int_t GetWindowTopX()
Returns current top x position of window on screen.
Definition TCanvas.cxx:1206
virtual void ToggleEventStatus()
Toggle event statusbar.
Definition TCanvas.cxx:2428
void Destructor()
Actual canvas destructor.
Definition TCanvas.cxx:692
Bool_t fUpdated
! Set to True when Update method was called
Definition TCanvas.h:65
void DeleteCanvasPainter()
assert on IsBatch() == false?
Definition TCanvas.cxx:2648
TPad * fPadSave
! Pointer to saved pad in HandleInput
Definition TCanvas.h:57
static Bool_t SupportAlpha()
Static function returning "true" if transparency is supported.
Definition TCanvas.cxx:2473
Bool_t fBatch
! True when in batchmode
Definition TCanvas.h:60
Bool_t fUseGL
! True when rendering is with GL
Definition TCanvas.h:63
Int_t fEventX
! Last X mouse position in canvas
Definition TCanvas.h:47
Bool_t IsBatch() const override
Definition TCanvas.h:173
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:906
Size_t fXsizeReal
Current size of canvas along X in CM.
Definition TCanvas.h:36
Bool_t HasMenuBar() const
Definition TCanvas.h:170
TVirtualPadPainter * GetCanvasPainter()
Access and (probably) creation of pad painter.
Definition TCanvas.cxx:2615
virtual void HighlightConnect(const char *slot)
This is "simplification" for function TCanvas::Connect with Highlighted signal for specific slot.
Definition TCanvas.cxx:1626
TPad * Pick(Int_t px, Int_t py, TObjLink *&pickobj) override
Search for an object at pixel position px,py.
Definition TCanvas.h:185
void Close(Option_t *option="") override
Close canvas.
Definition TCanvas.cxx:785
void SetFixedAspectRatio(Bool_t fixed=kTRUE) override
Fix canvas aspect ratio to current value if fixed is true.
Definition TCanvas.cxx:1989
virtual void Resize(Option_t *option="")
Recompute canvas parameters following a X11 Resize.
Definition TCanvas.cxx:1663
Color_t GetHighLightColor() const override
Definition TCanvas.h:140
Bool_t GetShowToolBar() const
Definition TCanvas.h:151
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:977
Bool_t IsFolder() const override
Is folder ?
Definition TCanvas.cxx:1483
Bool_t EnsurePSPainter(Bool_t create, TVirtualPadPainter *&oldp)
Replace canvas painter For intenral use only - when creating PS images.
Definition TCanvas.cxx:2625
UInt_t fWindowWidth
Width of window (including borders, etc.)
Definition TCanvas.h:42
TVirtualPadPainter * fPainter
! Canvas (pad) painter.
Definition TCanvas.h:67
void CopyPixmaps() override
Copy the canvas pixmap of the pad to the canvas.
Definition TCanvas.cxx:835
Bool_t IsGrayscale()
Check whether this canvas is to be drawn in grayscale mode.
Definition TCanvas.cxx:2566
TPad * fClickSelectedPad
! Pad containing currently click-selected object
Definition TCanvas.h:56
Bool_t fUpdating
! True when Updating the canvas
Definition TCanvas.h:61
void SaveSource(const char *filename="", Option_t *option="")
Save primitives in this canvas as a C++ macro file.
Definition TCanvas.cxx:1809
void SetCanvasImp(TCanvasImp *i)
Set canvas implementation If web-based implementation provided, some internal fields also initialized...
Definition TCanvas.cxx:2158
Color_t fHighLightColor
Highlight color of active pad.
Definition TCanvas.h:38
virtual void Size(Float_t xsizeuser=0, Float_t ysizeuser=0)
Set the canvas scale in centimeters.
Definition TCanvas.cxx:2198
virtual void ProcessedEvent(Int_t event, Int_t x, Int_t y, TObject *selected)
Emit ProcessedEvent() signal.
Definition TCanvas.cxx:1648
virtual void HandleInput(EEventType button, Int_t x, Int_t y)
Handle Input Events.
Definition TCanvas.cxx:1230
void UpdateAsync() override
Asynchronous pad update.
Definition TCanvas.cxx:2544
Size_t fXsizeUser
User specified size of canvas along X in CM.
Definition TCanvas.h:34
Int_t fEventY
! Last Y mouse position in canvas
Definition TCanvas.h:48
TVirtualPad * cd(Int_t subpadnumber=0) override
Set current canvas & pad.
Definition TCanvas.cxx:716
UInt_t fWindowHeight
Height of window (including menubar, borders, etc.)
Definition TCanvas.h:43
Int_t GetWindowTopY()
Returns current top y position of window on screen.
Definition TCanvas.cxx:1217
TObject * fClickSelected
! Currently click-selected object
Definition TCanvas.h:51
void SetCanvasSize(UInt_t ww, UInt_t wh) override
Set Width and Height of canvas to ww and wh respectively.
Definition TCanvas.cxx:1947
void Show()
Show canvas.
Definition TCanvas.cxx:2209
TPad * fSelectedPad
! Pad containing currently selected object
Definition TCanvas.h:55
virtual void Selected(TVirtualPad *pad, TObject *obj, Int_t event)
Emit Selected() signal.
Definition TCanvas.cxx:1634
Int_t fSelectedX
! X of selected object
Definition TCanvas.h:52
virtual void ToggleEditor()
Toggle editor.
Definition TCanvas.cxx:2450
TVirtualPad * GetSelectedPad() const override
Definition TCanvas.h:148
virtual void Picked(TPad *selpad, TObject *selected, Int_t event)
Emit Picked() signal.
Definition TCanvas.cxx:1588
TObject * fSelected
! Currently selected object
Definition TCanvas.h:50
void SavePrimitive(std::ostream &out, Option_t *option="") override
Save primitives in this canvas in C++ macro file with GUI.
Definition TCanvas.cxx:1775
void SetCursor(ECursor cursor) override
Set cursor.
Definition TCanvas.cxx:1961
void Streamer(TBuffer &) override
Stream a class object.
Definition TCanvas.cxx:2218
Bool_t GetShowToolTips() const
Definition TCanvas.h:153
Int_t fCanvasID
! Canvas identifier
Definition TCanvas.h:49
void SetGrayscale(Bool_t set=kTRUE)
Set whether this canvas should be painted in grayscale, and re-paint it if necessary.
Definition TCanvas.cxx:2575
void SetTitle(const char *title="") override
Set canvas title.
Definition TCanvas.cxx:2128
UInt_t fCh
Height of the canvas along Y (pixels)
Definition TCanvas.h:45
TContextMenu * fContextMenu
! Context menu pointer
Definition TCanvas.h:59
TAttCanvas fCatt
Canvas attributes.
Definition TCanvas.h:32
void SetName(const char *name="") override
Set canvas name.
Definition TCanvas.cxx:2020
UInt_t GetWindowWidth() const
Definition TCanvas.h:163
Bool_t fRetained
Retain structure flag.
Definition TCanvas.h:62
void DisconnectWidget()
Used by friend class TCanvasImp.
Definition TCanvas.cxx:2557
void FeedbackMode(Bool_t set)
Turn rubberband feedback mode on or off.
Definition TCanvas.cxx:1127
void ls(Option_t *option="") const override
List all pads.
Definition TCanvas.cxx:1499
void RaiseWindow()
Raise canvas window.
Definition TCanvas.cxx:1741
void Build()
Build a canvas. Called by all constructors.
Definition TCanvas.cxx:586
static Bool_t SaveAll(const std::vector< TPad * > &={}, const char *filename="", Option_t *option="")
Save provided pads/canvases into the image file(s) Filename can include printf argument for image num...
Definition TCanvas.cxx:2676
Int_t fWindowTopY
Top Y position of window (in pixels)
Definition TCanvas.h:41
void Paint(Option_t *option="") override
Paint canvas.
Definition TCanvas.cxx:1538
@ kResizeOpaque
Definition TCanvas.h:97
@ kShowToolTips
Definition TCanvas.h:99
@ kShowToolBar
Definition TCanvas.h:94
@ kMoveOpaque
Definition TCanvas.h:96
@ kIsGrayscale
Definition TCanvas.h:98
@ kShowEventStatus
Definition TCanvas.h:91
@ kAutoExec
Definition TCanvas.h:92
@ kMenuBar
Definition TCanvas.h:93
@ kShowEditor
Definition TCanvas.h:95
TClass * IsA() const override
Definition TCanvas.h:240
void Update() override
Update canvas pad buffers.
Definition TCanvas.cxx:2488
void ExecuteEvent(Int_t event, Int_t px, Int_t py) override
Execute action corresponding to one event.
Definition TCanvas.cxx:1109
void RunAutoExec()
Execute the list of TExecs in the current pad.
Definition TCanvas.cxx:1764
virtual void Cleared(TVirtualPad *pad)
Emit pad Cleared signal.
Definition TCanvas.cxx:767
UInt_t GetWw() const override
Definition TCanvas.h:165
TCanvasImp * fCanvasImp
! Window system specific canvas implementation
Definition TCanvas.h:58
UInt_t GetWh() const override
Definition TCanvas.h:166
virtual void Highlighted(TVirtualPad *pad, TObject *obj, Int_t x, Int_t y)
Emit Highlighted() signal.
Definition TCanvas.cxx:1607
void Flush()
Flush canvas buffers.
Definition TCanvas.cxx:1141
Size_t fYsizeUser
User specified size of canvas along Y in CM.
Definition TCanvas.h:35
Int_t fDoubleBuffer
Double buffer flag (0=off, 1=on)
Definition TCanvas.h:39
void ForceUpdate()
Force canvas update.
Definition TCanvas.cxx:1173
void CreatePainter()
Probably, TPadPainter must be placed in a separate ROOT module - "padpainter" (the same as "histpaint...
Definition TCanvas.cxx:2594
void SetSelected(TObject *obj) override
Set selected canvas.
Definition TCanvas.cxx:2119
void MoveOpaque(Int_t set=1)
Set option to move objects/pads in a canvas.
Definition TCanvas.cxx:1530
static Bool_t fgIsFolder
Indicates if canvas can be browsed as a folder.
Definition TCanvas.h:69
void Closed() override
Emit Closed signal.
Definition TCanvas.cxx:775
Bool_t IsWeb() const override
Is web canvas.
Definition TCanvas.cxx:1491
void SetWindowPosition(Int_t x, Int_t y)
Set canvas window position.
Definition TCanvas.cxx:2137
TString fDISPLAY
Name of destination screen.
Definition TCanvas.h:33
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:2052
Int_t fEvent
! Type of current or last handled event
Definition TCanvas.h:46
Bool_t GetShowEventStatus() const
Definition TCanvas.h:150
TString fSelectedOpt
! Drawing option of selected object
Definition TCanvas.h:54
Int_t fSelectedY
! Y of selected object
Definition TCanvas.h:53
Bool_t fDrawn
! Set to True when the Draw method is called
Definition TCanvas.h:64
void SetBatch(Bool_t batch=kTRUE) override
Toggle batch mode.
Definition TCanvas.cxx:1929
Bool_t UseGL() const
Definition TCanvas.h:230
void ResizeOpaque(Int_t set=1)
Set option to resize objects/pads in a canvas.
Definition TCanvas.cxx:1756
virtual void ToggleToolBar()
Toggle toolbar.
Definition TCanvas.cxx:2439
virtual TObject * DrawClonePad()
Draw a clone of this canvas into the current pad In an interactive session, select the destination/cu...
Definition TCanvas.cxx:923
@ kRealNew
Definition TClass.h:110
static ENewType IsCallingNew()
Static method returning the defConstructor flag passed to TClass::New().
Definition TClass.cxx:5960
void Browse(TBrowser *b) override
Browse this collection (called by TBrowser).
The color creation and management class.
Definition TColor.h:22
static const TArrayI & GetPalette()
Static function returning the current active palette.
Definition TColor.cxx:1521
static void SaveColorsPalette(std::ostream &out)
Store current palette in the output macro.
Definition TColor.cxx:3649
static TClass * Class()
static TString SavePrimitiveColor(Int_t ci)
Convert color in C++ statement which can be used in SetColor directives Produced statement either inc...
Definition TColor.cxx:2556
static Bool_t DefinedColors(Int_t set_always_on=0)
Static method returning kTRUE if some new colors have been defined after initialisation or since the ...
Definition TColor.cxx:1542
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 * AsString() const
Return the date & time as a string (ctime() format).
Definition TDatime.cxx:101
static TClass * Class()
virtual Int_t GetValue(const char *name, Int_t dflt) const
Returns the integer value for a resource.
Definition TEnv.cxx:503
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:109
Option_t * GetOption() const
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:708
void Add(TObject *obj) override
Definition TList.h:81
TObject * Remove(TObject *obj) override
Remove object from the list.
Definition TList.cxx:952
An array of TObjects.
Definition TObjArray.h:31
static TClass * Class()
Mother of all ROOT objects.
Definition TObject.h:42
virtual void Clear(Option_t *="")
Definition TObject.h:127
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition TObject.h:204
virtual TObject * Clone(const char *newname="") const
Make a clone of an object using the Streamer facility.
Definition TObject.cxx:240
R__ALWAYS_INLINE Bool_t IsOnHeap() const
Definition TObject.h:160
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:224
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:1081
virtual void ExecuteEvent(Int_t event, Int_t px, Int_t py)
Execute action corresponding to an event at (px,py).
Definition TObject.cxx:412
virtual void Delete(Option_t *option="")
Delete this object.
Definition TObject.cxx:265
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:885
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition TObject.cxx:546
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:1095
virtual void Pop()
Pop on object drawn in a pad to the top of the display list.
Definition TObject.cxx:637
@ kNoContextMenu
if object does not want context menu
Definition TObject.h:78
@ kMustCleanup
if object destructor must call RecursiveRemove()
Definition TObject.h:73
virtual void Info(const char *method, const char *msgfmt,...) const
Issue info message.
Definition TObject.cxx:1069
static TClass * Class()
Implement TVirtualPadPainter which abstracts painting operations.
Definition TPadPainter.h:27
The most important graphics class in the ROOT system.
Definition TPad.h:28
Short_t GetBorderMode() const override
Definition TPad.h:200
void SetBorderSize(Short_t bordersize) override
Definition TPad.h:331
Int_t GetTicky() const override
Definition TPad.h:240
void PaintBorder(Color_t color, Bool_t tops)
Paint the pad border.
Definition TPad.cxx:3762
void ResizePad(Option_t *option="") override
Compute pad conversion coefficients.
Definition TPad.cxx:5752
void SetGrid(Int_t valuex=1, Int_t valuey=1) override
Definition TPad.h:340
void SavePrimitive(std::ostream &out, Option_t *option="") override
Save primitives in this pad on the C++ source file out.
Definition TPad.cxx:5969
Bool_t GetGridx() const override
Definition TPad.h:236
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:6182
Double_t GetY2() const override
Definition TPad.h:244
void Close(Option_t *option="") override
Delete all primitives in pad and pad itself.
Definition TPad.cxx:1087
void PaintModified() override
Traverse pad hierarchy and (re)paint only modified pads.
Definition TPad.cxx:3878
const char * GetTitle() const override
Returns title of object.
Definition TPad.h:262
Double_t fX1
X of lower X coordinate.
Definition TPad.h:36
TList * GetListOfPrimitives() const override
Definition TPad.h:246
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:6242
void UseCurrentStyle() override
Force a copy of current style for all objects in pad.
Definition TPad.cxx:7011
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:5459
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:722
Int_t GetTickx() const override
Definition TPad.h:239
Double_t fAspectRatio
ratio of w/h in case of fixed ratio
Definition TPad.h:82
void SetLogy(Int_t value=1) override
Set Lin/Log scale for Y.
Definition TPad.cxx:6171
TCanvas * fCanvas
! Pointer to mother canvas
Definition TPad.h:106
Bool_t fFixedAspectRatio
True if fixed aspect ratio.
Definition TPad.h:104
void Modified(Bool_t flag=true) override
Mark pad modified Will be repainted when TCanvas::Update() will be called next time.
Definition TPad.cxx:7472
void ls(Option_t *option="") const override
List all primitives in pad.
Definition TPad.cxx:3198
void Streamer(TBuffer &) override
Stream a class object.
Definition TPad.cxx:6797
TString fTitle
Pad title.
Definition TPad.h:110
void CopyPixmaps() override
Copy the sub-pixmaps of the pad to the canvas.
Definition TPad.cxx:1159
void CopyPixmap() override
Copy the pixmap of the pad to the canvas.
Definition TPad.cxx:1146
Double_t GetY1() const override
Definition TPad.h:243
Bool_t GetGridy() const override
Definition TPad.h:237
void ExecuteEvent(Int_t event, Int_t px, Int_t py) override
Execute action corresponding to one event.
Definition TPad.cxx:1951
Int_t GetLogz() const override
Definition TPad.h:259
Short_t GetBorderSize() const override
Definition TPad.h:201
TList * fPrimitives
->List of primitives (subpads)
Definition TPad.h:107
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:3692
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:2807
TVirtualPad * cd(Int_t subpadnumber=0) override
Set Current pad.
Definition TPad.cxx:693
Int_t GetLogy() const override
Definition TPad.h:258
void SetBorderMode(Short_t bordermode) override
Definition TPad.h:330
void SetTicks(Int_t valuex=1, Int_t valuey=1) override
Definition TPad.h:360
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:6157
Int_t GetLogx() const override
Definition TPad.h:257
Double_t GetX2() const override
Definition TPad.h:242
Double_t GetX1() const override
Definition TPad.h:241
TPad * fMother
! pointer to mother of the list
Definition TPad.h:105
const char * GetName() const override
Returns name of object.
Definition TPad.h:261
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:865
static Int_t IncreaseDirLevel()
Increase the indentation level for ls().
Definition TROOT.cxx:2897
static void IndentLevel()
Functions used by ls() to indent an object hierarchy.
Definition TROOT.cxx:2905
static Int_t DecreaseDirLevel()
Decrease the indentation level for ls().
Definition TROOT.cxx:2754
Basic string class.
Definition TString.h:138
void ToLower()
Change string to lower-case.
Definition TString.cxx:1189
TString & ReplaceSpecialCppChars()
Find special characters which are typically used in printf() calls and replace them by appropriate es...
Definition TString.cxx:1121
Bool_t IsDigit() const
Returns true if all characters in string are digits (0-9) or white spaces, i.e.
Definition TString.cxx:1836
@ kBoth
Definition TString.h:284
virtual void Streamer(TBuffer &)
Stream a string object.
Definition TString.cxx:1418
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString.
Definition TString.cxx:2384
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:641
Double_t GetTimeOffset() const
Definition TStyle.h:271
Int_t GetOptLogy() const
Definition TStyle.h:250
Int_t GetOptStat() const
Definition TStyle.h:247
Int_t GetOptTitle() const
Definition TStyle.h:248
void SetCanvasBorderSize(Width_t size=1)
Definition TStyle.h:348
Float_t GetScreenFactor() const
Definition TStyle.h:258
Int_t GetPadTickX() const
Definition TStyle.h:219
Bool_t IsReading() const
Definition TStyle.h:300
void SetCanvasColor(Color_t color=19)
Definition TStyle.h:347
Float_t GetPadRightMargin() const
Definition TStyle.h:216
void SetCanvasBorderMode(Int_t mode=1)
Definition TStyle.h:349
Int_t GetCanvasDefH() const
Definition TStyle.h:193
Int_t GetCanvasDefX() const
Definition TStyle.h:195
Bool_t GetPadGridY() const
Definition TStyle.h:218
Float_t GetPadLeftMargin() const
Definition TStyle.h:215
void SetPalette(Int_t ncolors=kBird, Int_t *colors=nullptr, Float_t alpha=1.)
See TColor::SetPalette.
Definition TStyle.cxx:1889
Bool_t GetCanvasPreferGL() const
Definition TStyle.h:189
Int_t GetCanvasDefY() const
Definition TStyle.h:196
Bool_t GetPadGridX() const
Definition TStyle.h:217
Int_t GetPadTickY() const
Definition TStyle.h:220
Color_t GetCanvasColor() const
Definition TStyle.h:190
Float_t GetPadBottomMargin() const
Definition TStyle.h:213
Int_t GetCanvasDefW() const
Definition TStyle.h:194
Int_t GetOptLogx() const
Definition TStyle.h:249
Int_t GetCanvasBorderMode() const
Definition TStyle.h:192
Width_t GetCanvasBorderSize() const
Definition TStyle.h:191
Int_t GetOptFit() const
Definition TStyle.h:246
Int_t GetOptLogz() const
Definition TStyle.h:251
Float_t GetPadTopMargin() const
Definition TStyle.h:214
virtual Bool_t ExpandPathName(TString &path)
Expand a pathname getting rid of special shell characters like ~.
Definition TSystem.cxx:1285
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 SetAttFill(const TAttFill &att)
Set fill attributes.
virtual void ClearDrawable()=0
virtual void UpdateDrawable(Int_t)
virtual void SetAttMarker(const TAttMarker &att)
Set marker attributes.
virtual void LockPainter()
Empty definition.
virtual void SetDrawMode(Int_t, Int_t)
virtual void SetCursor(Int_t win, ECursor cursor)
Set cursor for specified device, redirect to gVirtualX.
virtual TClass * IsA() const
virtual void SetAttLine(const TAttLine &att)
Set line attributes.
static TVirtualPadPainter * PadPainter(Option_t *opt="")
Create a pad painter of specified type.
virtual void SetAttText(const TAttText &att)
Set text attributes.
virtual void SelectDrawable(Int_t device)=0
virtual void InitPainter()
Empty definition.
virtual void SetDoubleBuffer(Int_t device, Int_t mode)
Set double buffer mode for specified device, redirect to gVirtualX.
small helper class to store/restore gPad context in TPad methods
Definition TVirtualPad.h:61
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:409
Int_t Nint(T x)
Round to nearest integer. Rounds half integers to the nearest even integer.
Definition TMath.h:704
th1 Draw()