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