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