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