Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TStyle.cxx
Go to the documentation of this file.
1// @(#)root/base:$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 <cstdio>
14#include <cctype>
15#include <cmath>
16#include <iostream>
17#include <fstream>
18
19#include "strlcpy.h"
20#include "TApplication.h"
21#include "TColor.h"
22#include "TDatime.h"
23#include "TROOT.h"
24#include "TStyle.h"
25#include "TSystem.h"
26#include "TVirtualPad.h"
27#include "TVirtualMutex.h"
28#include "TEnv.h"
29
31const UInt_t kTakeStyle = BIT(17);
32
33
34/** \class TStyle
35\ingroup Base
36 \ingroup GraphicsAtt
37
38TStyle objects may be created to define special styles.
39By default ROOT creates a default style that can be accessed via
40the gStyle pointer.
41
42This class includes functions to set some of the following object attributes.
43 - Canvas
44 - Pad
45 - Histogram axis
46 - Lines
47 - Fill areas
48 - Text
49 - Markers
50 - Functions
51 - Histogram Statistics and Titles
52
53All objects that can be drawn in a pad inherit from one or more attribute classes
54like TAttLine, TAttFill, TAttText, TAttMarker. When the objects are created, their
55default attributes are taken from the current style. The current style is an object
56of the class[TStyle](https://root.cern/doc/master/classTStyle.html) and can be
57referenced via the global variable `gStyle` (in TStyle.h).
58
59ROOT provides two styles called "Default" and "Plain". The "Default"
60style is created simply by:
61
62~~~ {.cpp}
63 auto default = new TStyle("Default","Default Style");
64~~~
65
66The `Plain` style can be used if you are working on a monochrome display or
67if you want to get a "conventional" Postscript output. These are the instructions
68in the ROOT constructor to create the `Plain` style.
69
70```
71auto plain = new TStyle("Plain","Plain Style (no colors/fill areas)");
72
73 plain->SetCanvasBorderMode(0);
74 plain->SetPadBorderMode(0);
75 plain->SetPadColor(0);
76 plain->SetCanvasColor(0);
77 plain->SetTitleColor(0);
78 plain->SetStatColor(0);
79```
80
81You can set the current style with:
82
83```
84gROOT->SetStyle(style_name);
85```
86
87You can get a pointer to an existing style with:
88
89```
90auto style = gROOT->GetStyle(style_name);
91```
92
93You can create additional styles with:
94
95```
96 TStyle *st1 = new TStyle("st1","my style");
97 st1->Set....
98 st1->cd(); this becomes now the current style gStyle
99```
100
101In your [rootlogon.C](https://root.cern/doc/master/classexamples/startsession.log.html)
102file, you can redefine the default parameters via statements like:
103
104```
105 gStyle->SetStatX(0.7);
106 gStyle->SetStatW(0.2);
107 gStyle->SetLabelOffset(1.2);
108 gStyle->SetLabelFont(72);
109```
110
111Note that when an object is created, its attributes are taken from the current
112style. For example, you may have created an histogram in a previous session,
113saved it in a file. Meanwhile, if you have changed the style, the histogram will
114be drawn with the old attributes. You can force the current style attributes to
115be set when you read an object from a file by calling:
116
117```
118gROOT->ForceStyle();
119```
120
121before reading the objects from the file.
122
123Let's assume you have a canvas or pad with your histogram or any other object,
124you can force these objects to get the attributes of the current style via:
125
126```
127canvas->UseCurrentStyle();
128```
129
130The description of the style functions should be clear from the name of the
131TStyle Setters or Getters. Some functions have an extended description, in particular:
132
133 - TStyle:SetLabelFont.
134 - TStyle:SetLineStyleString, to set the format of dashed lines.
135 - TStyle:SetOptStat.
136 - TStyle:SetPalette to change the colors palette.
137 - TStyle:SetTitleOffset.
138
139*/
140
141////////////////////////////////////////////////////////////////////////////////
142/// Default constructor.
143
145{
146 Reset();
147}
148
149////////////////////////////////////////////////////////////////////////////////
150/// Create a new TStyle.
151/// The following names are reserved to create special styles:
152///
153/// - `Classic`: Similar to `Default` style set in TStyle::Reset
154/// - `Plain`: a black&white oriented style
155/// - `Bold`
156/// - `Video`
157/// - `Pub`
158/// - `Modern`: Used when ROOT starts
159/// - `ATLAS`: style used by the ATLAS experiment
160/// - `BELLE2`: style used by the BELLE II experiment
161/// (see the definition of these styles below).
162///
163/// Note a side-effect of calling gStyle->SetFillColor(0). This is nearly
164/// equivalent of selecting the "Plain" style.
165///
166/// Many graphics attributes may be set via the TStyle, see in particular
167/// - TStyle::SetNdivisions
168/// - TStyle::SetAxisColor
169/// - TStyle::SetHeaderPS
170/// - TStyle::SetTitlePS
171/// - TStyle::SetLabelColor
172/// - TStyle::SetLabelFont
173/// - TStyle::SetLabelOffset
174/// - TStyle::SetLabelSize
175/// - TStyle::SetOptDate
176/// - TStyle::SetLineStyleString
177/// - TStyle::SetOptFit
178/// - TStyle::SetOptStat
179/// - TStyle::SetPaperSize
180/// - TStyle::SetTickLength
181/// - TStyle::SetTitleOffset
182/// - TStyle::SetTitleSize
183/// - TStyle::SetPalette
184/// - TStyle::SetTimeOffset
185/// - TStyle::SetStripDecimals
186///
187/// The current style is pointed by gStyle.
188///
189/// When calling myStyle->cd(), gStyle is set to myStyle.
190///
191/// One can also use gROOT to change the current style, e.g.
192///
193/// gROOT->SetStyle("Plain") will change the current style gStyle to the
194/// "Plain" style
195///
196/// See also TROOT::ForceStyle and TROOT::UseCurrentStyle
197
198TStyle::TStyle(const char *name, const char *title)
199{
201
202 SetNameTitle(style_name, title);
203
204 // If another style was already created with the same name, it is overwrite.
205 delete gROOT->GetStyle(style_name);
206
207 Reset();
208
209 {
211 gROOT->GetListOfStyles()->Add(this);
212 }
213
214 if (strcmp(style_name,"Modern") == 0) {
215 // Modern style
221 SetPadColor(0);
222 SetStatColor(0);
223 SetTitleFont(42,"");
224 SetLabelFont(42,"x");
225 SetTitleFont(42,"x");
226 SetLabelFont(42,"y");
227 SetTitleFont(42,"y");
228 SetLabelFont(42,"z");
229 SetTitleFont(42,"z");
230 SetStatFont(42);
231 SetLabelSize(0.035,"x");
232 SetTitleSize(0.035,"x");
233 SetLabelSize(0.035,"y");
234 SetTitleSize(0.035,"y");
235 SetLabelSize(0.035,"z");
236 SetTitleSize(0.035,"z");
237 SetTitleSize(0.050,"");
238 SetTitleAlign(23);
239 SetTitleX(0.5);
242 SetTitleStyle(0);
243 SetTitleOffset(0.,"Y");
245 SetOptStat(1111);
246 SetStatY(0.935);
250 SetLegendFillStyle(1001);
251 SetLegendFont(42);
253 SetFuncWidth(2);
254 SetFuncColor(2);
255 }
256 if (strcmp(style_name,"Plain") == 0) {
257 // May be a standard style to be initialised
262 SetPadColor(0);
266 SetStatColor(0);
269 return;
270 }
271 if (strcmp(style_name,"Bold") == 0) {
272 // Authors: Art Poskanzer and Jim Thomas, LBNL, Oct. 2000
273 SetPalette(1,nullptr);
274 SetCanvasColor(10);
278 SetPadColor(10);
279 SetPadTickX(1);
280 SetPadTickY(1);
281 SetPadBottomMargin(0.15);
282 SetPadLeftMargin(0.15);
285 SetFuncWidth(3);
287 SetLineWidth(3);
288 SetLabelSize(0.05,"xyz");
289 SetLabelOffset(0.01,"y");
290 SetLabelColor(kBlue,"xy");
291 SetTitleSize(0.06,"xyz");
292 SetTitleOffset(1.3,"Y");
295 SetStatColor(10);
296 return;
297 }
298 if (strcmp(style_name,"Video") == 0) {
299 // Author: Art Poskanzer, LBNL, Oct. 1999
300 SetPalette(1,nullptr);
301 SetCanvasColor(10);
305 SetPadColor(10);
306 SetPadTickX(1);
307 SetPadTickY(1);
309 SetPadLeftMargin(0.2);
312 SetLabelSize(0.06,"xyz");
313 SetLabelColor(kBlue,"xyz");
314 SetTitleSize(0.08,"xyz");
317 SetStatColor(10);
318 SetFuncWidth(8);
320 SetLineWidth(3);
321 return;
322 }
323 if (strcmp(style_name,"Pub") == 0) {
324 // Authors: Art Poskanzer and Jim Thomas, LBNL, Oct. 2000
325 SetOptTitle(0);
326 SetOptStat(0);
327 SetPalette(8,nullptr);
328 SetCanvasColor(10);
332 SetPadColor(10);
333 SetPadTickX(1);
334 SetPadTickY(1);
335 SetPadBottomMargin(0.15);
336 SetPadLeftMargin(0.15);
339 SetFuncWidth(3);
341 SetLineWidth(3);
342 SetLabelSize(0.05,"xyz");
343 SetLabelOffset(0.01,"y");
344 SetLabelColor(kBlack,"xyz");
345 SetTitleSize(0.06,"xyz");
346 SetTitleOffset(1.3,"y");
349 return;
350 }
351 if (strcmp(style_name,"ATLAS") == 0) {
352 // Author: M.Sutton - Atlas Collaboration 2010
358 SetPadColor(0);
359 SetStatColor(0);
360 SetPaperSize(20,26);
361 SetPadTopMargin(0.05);
362 SetPadRightMargin(0.05);
363 SetPadBottomMargin(0.16);
364 SetPadLeftMargin(0.16);
365 SetTitleXOffset(1.4);
366 SetTitleYOffset(1.4);
367 Int_t font = 42;
368 Double_t tsize=0.05;
369 SetTextFont(font);
371 SetLabelFont(font,"x");
372 SetTitleFont(font,"x");
373 SetLabelFont(font,"y");
374 SetTitleFont(font,"y");
375 SetLabelFont(font,"z");
376 SetTitleFont(font,"z");
377 SetLabelSize(tsize,"x");
378 SetTitleSize(tsize,"x");
379 SetLabelSize(tsize,"y");
380 SetTitleSize(tsize,"y");
381 SetLabelSize(tsize,"z");
382 SetTitleSize(tsize,"z");
383 SetMarkerStyle(20);
384 SetMarkerSize(1.2);
386 SetLineStyleString(2,"[12 12]");
387 SetEndErrorSize(0.); // get rid of error bar caps
388 SetOptTitle(0);
389 SetOptStat(0);
390 SetOptFit(0);
391 SetPadTickX(1);
392 SetPadTickY(1);
393 }
394 if (strcmp(style_name,"BELLE2") == 0) {
395 // use plain black on white colors
396 Int_t icol=0; // WHITE
404 //SetFillColor(icol); // don't use: white fill color for *all* objects
405
406 // set the paper & margin sizes
407 SetPaperSize(20,26);
408
409 // set margin sizes
410 SetPadTopMargin(0.05);
411 SetPadRightMargin(0.05);
412 SetPadBottomMargin(0.16);
413 SetPadLeftMargin(0.16);
414
415 // set title offsets (for axis label)
416 SetTitleXOffset(1.0);
417 SetTitleYOffset(1.0);
418
419 // use large fonts
420 //Int_t font=72; // Helvetica italics
421 Int_t font=42; // Helvetica
422 Double_t tsize=0.05;
423 SetTextFont(font);
425
426 SetLabelFont(font,"x");
427 SetTitleFont(font,"x");
428 SetLabelFont(font,"y");
429 SetTitleFont(font,"y");
430 SetLabelFont(font,"z");
431 SetTitleFont(font,"z");
432
433 SetLabelSize(tsize,"x");
434 SetTitleSize(.065,"x");
435 SetLabelSize(tsize,"y");
436 SetTitleSize(.065,"y");
437 SetLabelSize(tsize,"z");
438 SetTitleSize(.065,"z");
439
440 SetTitleOffset(1.1,"x");
441 SetTitleOffset(1.1,"y");
442 SetTitleOffset(1.1,"z");
443
444 SetLabelOffset(0.015,"x");
445 SetLabelOffset(0.015,"y");
446 SetLabelOffset(0.015,"z");
447
448 SetTickLength(0.03,"x");
449 SetTickLength(0.02,"y"); // This way we slightly achieve equal length ticks for x and y
450
451 // use bold lines and markers
452 SetMarkerStyle(20);
453 SetMarkerSize(0.9);
455 SetLineStyleString(2,"[12 12]"); // postscript dashes
456
457 // get rid of X error bars
458 SetErrorX(0.001);
459 // get rid of error bar caps
460 SetEndErrorSize(0.);
461
462 // do not display any of the standard histogram decorations
463 SetOptTitle(0);
464 SetOptStat(0);
465 SetOptFit(0);
466
467 // put tick marks on top and RHS of plots
468 SetPadTickX(0);
469 SetPadTickY(0);
470
472 }
473}
474
475////////////////////////////////////////////////////////////////////////////////
476/// Destructor.
477
479{
481 gROOT->GetListOfStyles()->Remove(this);
482 if (gStyle == this) gStyle = (TStyle*)gROOT->GetListOfStyles()->Last();
483}
484
485////////////////////////////////////////////////////////////////////////////////
486/// Copy constructor
487
489{
490 style.TStyle::Copy(*this);
491}
492
493////////////////////////////////////////////////////////////////////////////////
494/// Assignment operator.
495
497{
498 if (this != &style)
499 style.TStyle::Copy(*this);
500 return *this;
501}
502
503////////////////////////////////////////////////////////////////////////////////
504// Return axis number (1 for X, 2 for Y, 3 for Z), otherwise 0
505
507{
508 UChar_t a = axis ? *axis : 0;
509 a -= (a >= 'x') ? 'x' : 'X'; // toupper and a-='X'; intentional underflow
510 return (a > 2) ? 0 : (Int_t)(a+1);
511}
512
513////////////////////////////////////////////////////////////////////////////////
514/// Browse the style object.
515
517{
518 cd();
519}
520
521////////////////////////////////////////////////////////////////////////////////
522/// Create some standard styles.
523
525{
526 TColor *col = new TColor(); // force the initialisation of fgPalette
527 new TStyle("Plain", "Plain Style (no colors/fill areas)");
528 new TStyle("Bold", "Bold Style");
529 new TStyle("Video", "Style for video presentation histograms");
530 new TStyle("Pub", "Style for Publications");
531 new TStyle("Classic","Classic Style");
532 new TStyle("Default","Equivalent to Classic");
533 new TStyle("Modern", "Modern Style");
534 new TStyle("ATLAS", "ATLAS Style");
535 new TStyle("BELLE2", "Belle II Style");
536 delete col;
537}
538
539////////////////////////////////////////////////////////////////////////////////
540/// Change current style.
541
543{
544 gStyle = this;
545}
546
547////////////////////////////////////////////////////////////////////////////////
548/// Copy this style.
549
550void TStyle::Copy(TObject &obj) const
551{
552 TAttLine::Copy(((TStyle&)obj));
553 TAttFill::Copy(((TStyle&)obj));
554 TAttMarker::Copy(((TStyle&)obj));
555 TAttText::Copy(((TStyle&)obj));
556 fXaxis.Copy(((TStyle&)obj).fXaxis);
557 fYaxis.Copy(((TStyle&)obj).fYaxis);
558 fZaxis.Copy(((TStyle&)obj).fZaxis);
559 fAttDate.Copy(((TStyle&)obj).fAttDate);
560 ((TStyle&)obj).fIsReading = fIsReading;
561 ((TStyle&)obj).fScreenFactor = fScreenFactor;
562 ((TStyle&)obj).fCanvasPreferGL = fCanvasPreferGL;
563 ((TStyle&)obj).fCanvasColor = fCanvasColor;
564 ((TStyle&)obj).fCanvasBorderSize = fCanvasBorderSize;
565 ((TStyle&)obj).fCanvasBorderMode = fCanvasBorderMode;
566 ((TStyle&)obj).fCanvasDefH = fCanvasDefH;
567 ((TStyle&)obj).fCanvasDefW = fCanvasDefW;
568 ((TStyle&)obj).fCanvasDefX = fCanvasDefX;
569 ((TStyle&)obj).fCanvasDefY = fCanvasDefY;
570 ((TStyle&)obj).fPadColor = fPadColor;
571 ((TStyle&)obj).fPadBorderSize = fPadBorderSize;
572 ((TStyle&)obj).fPadBorderMode = fPadBorderMode;
573 ((TStyle&)obj).fPadBottomMargin = fPadBottomMargin;
574 ((TStyle&)obj).fPadTopMargin = fPadTopMargin;
575 ((TStyle&)obj).fPadLeftMargin = fPadLeftMargin;
576 ((TStyle&)obj).fPadRightMargin = fPadRightMargin;
577 ((TStyle&)obj).fPadGridX = fPadGridX;
578 ((TStyle&)obj).fPadGridY = fPadGridY;
579 ((TStyle&)obj).fPadTickX = fPadTickX;
580 ((TStyle&)obj).fPadTickY = fPadTickY;
581 ((TStyle&)obj).fPaperSizeX = fPaperSizeX;
582 ((TStyle&)obj).fPaperSizeY = fPaperSizeY;
583 ((TStyle&)obj).fFuncColor = fFuncColor;
584 ((TStyle&)obj).fFuncStyle = fFuncStyle;
585 ((TStyle&)obj).fFuncWidth = fFuncWidth;
586 ((TStyle&)obj).fGridColor = fGridColor;
587 ((TStyle&)obj).fGridStyle = fGridStyle;
588 ((TStyle&)obj).fGridWidth = fGridWidth;
589 ((TStyle&)obj).fHatchesSpacing = fHatchesSpacing;
590 ((TStyle&)obj).fHatchesLineWidth = fHatchesLineWidth;
591 ((TStyle&)obj).fFrameFillColor = fFrameFillColor;
592 ((TStyle&)obj).fFrameFillStyle = fFrameFillStyle;
593 ((TStyle&)obj).fFrameLineColor = fFrameLineColor;
594 ((TStyle&)obj).fFrameLineStyle = fFrameLineStyle;
595 ((TStyle&)obj).fFrameLineWidth = fFrameLineWidth;
596 ((TStyle&)obj).fFrameBorderSize = fFrameBorderSize;
597 ((TStyle&)obj).fFrameBorderMode = fFrameBorderMode;
598 ((TStyle&)obj).fHistFillColor = fHistFillColor;
599 ((TStyle&)obj).fHistFillStyle = fHistFillStyle;
600 ((TStyle&)obj).fHistLineColor = fHistLineColor;
601 ((TStyle&)obj).fHistLineStyle = fHistLineStyle;
602 ((TStyle&)obj).fHistLineWidth = fHistLineWidth;
603 ((TStyle&)obj).fHistMinimumZero = fHistMinimumZero;
604 ((TStyle&)obj).fHistTopMargin = fHistTopMargin;
605 ((TStyle&)obj).fBarWidth = fBarWidth;
606 ((TStyle&)obj).fBarOffset = fBarOffset;
607 ((TStyle&)obj).fDrawBorder = fDrawBorder;
608 ((TStyle&)obj).fOptLogx = fOptLogx;
609 ((TStyle&)obj).fOptLogy = fOptLogy;
610 ((TStyle&)obj).fOptLogz = fOptLogz;
611 ((TStyle&)obj).fOptDate = fOptDate;
612 ((TStyle&)obj).fOptFile = fOptFile;
613 ((TStyle&)obj).fOptFit = fOptFit;
614 ((TStyle&)obj).fOptStat = fOptStat;
615 ((TStyle&)obj).fOptTitle = fOptTitle;
616 ((TStyle&)obj).fEndErrorSize = fEndErrorSize;
617 ((TStyle&)obj).fErrorX = fErrorX;
618 ((TStyle&)obj).fStatColor = fStatColor;
619 ((TStyle&)obj).fStatTextColor = fStatTextColor;
620 ((TStyle&)obj).fStatBorderSize = fStatBorderSize;
621 ((TStyle&)obj).fStatFont = fStatFont;
622 ((TStyle&)obj).fStatFontSize = fStatFontSize;
623 ((TStyle&)obj).fStatStyle = fStatStyle;
624 ((TStyle&)obj).fStatFormat = fStatFormat;
625 ((TStyle&)obj).fStatW = fStatW;
626 ((TStyle&)obj).fStatH = fStatH ;
627 ((TStyle&)obj).fStatX = fStatX;
628 ((TStyle&)obj).fStatY = fStatY;
629 ((TStyle&)obj).fTitleAlign = fTitleAlign;
630 ((TStyle&)obj).fTitleColor = fTitleColor;
631 ((TStyle&)obj).fTitleTextColor = fTitleTextColor;
632 ((TStyle&)obj).fTitleFont = fTitleFont;
633 ((TStyle&)obj).fTitleFontSize = fTitleFontSize;
634 ((TStyle&)obj).fTitleStyle = fTitleStyle;
635 ((TStyle&)obj).fTitleBorderSize = fTitleBorderSize;
636 ((TStyle&)obj).fTitleW = fTitleW;
637 ((TStyle&)obj).fTitleH = fTitleH;
638 ((TStyle&)obj).fTitleX = fTitleX;
639 ((TStyle&)obj).fTitleY = fTitleY;
640 ((TStyle&)obj).fDateX = fDateX;
641 ((TStyle&)obj).fDateY = fDateY;
642 ((TStyle&)obj).fFitFormat = fFitFormat;
643 ((TStyle&)obj).fPaintTextFormat = fPaintTextFormat;
644 ((TStyle&)obj).fShowEventStatus = fShowEventStatus;
645 ((TStyle&)obj).fShowEditor = fShowEditor;
646 ((TStyle&)obj).fShowToolBar = fShowToolBar;
647 ((TStyle&)obj).fLegoInnerR = fLegoInnerR;
648 ((TStyle&)obj).fStripDecimals = fStripDecimals;
649 ((TStyle&)obj).fNumberContours = fNumberContours;
650 ((TStyle&)obj).fLegendBorderSize = fLegendBorderSize;
651 ((TStyle&)obj).fLegendFillColor = fLegendFillColor;
652 ((TStyle&)obj).fLegendFillStyle = fLegendFillStyle;
653 ((TStyle&)obj).fLegendFont = fLegendFont;
654 ((TStyle&)obj).fLegendTextSize = fLegendTextSize;
655
656 for (Int_t i=0;i<30;i++)
657 ((TStyle&)obj).fLineStyle[i] = fLineStyle[i];
658
659 ((TStyle&)obj).fHeaderPS = fHeaderPS;
660 ((TStyle&)obj).fTitlePS = fTitlePS;
661 ((TStyle&)obj).fLineScalePS = fLineScalePS;
662 ((TStyle&)obj).fJoinLinePS = fJoinLinePS;
663 ((TStyle&)obj).fCapLinePS = fCapLinePS;
664 ((TStyle&)obj).fColorModelPS = fColorModelPS;
665 ((TStyle&)obj).fTimeOffset = fTimeOffset;
666 ((TStyle&)obj).fImageScaling = fImageScaling;
667
668 ((TStyle&)obj).fCandleWhiskerRange = fCandleWhiskerRange;
669 ((TStyle&)obj).fCandleBoxRange = fCandleBoxRange;
670 ((TStyle&)obj).fCandleScaled = fCandleScaled;
671 ((TStyle&)obj).fViolinScaled = fViolinScaled;
672
673 ((TStyle&)obj).fOrthoCamera = fOrthoCamera;
674
675 ((TStyle&)obj).fXAxisExpXOffset = fXAxisExpXOffset;
676 ((TStyle&)obj).fXAxisExpYOffset = fXAxisExpYOffset;
677 ((TStyle&)obj).fYAxisExpXOffset = fYAxisExpXOffset;
678 ((TStyle&)obj).fYAxisExpYOffset = fYAxisExpYOffset;
679 ((TStyle&)obj).fAxisMaxDigits = fAxisMaxDigits;
680}
681
682////////////////////////////////////////////////////////////////////////////////
683/// Function used by the TStyle manager when drawing a canvas showing the
684/// current style.
685
687{
688 gPad->SetSelected(this);
689 return 0;
690}
691
692////////////////////////////////////////////////////////////////////////////////
693/// Reset.
694
696{
702 SetFillStyle(1001);
703 SetFillColor(19);
704 fXaxis.ResetAttAxis("X");
705 fYaxis.ResetAttAxis("Y");
706 fZaxis.ResetAttAxis("Z");
707 if (gEnv) fCanvasPreferGL = gEnv->GetValue("OpenGL.CanvasPreferGL",0);
708 else fCanvasPreferGL = kFALSE;
709 fCanvasColor = 19;
712 fCanvasDefH = 500;
713 fCanvasDefW = 700;
714 fCanvasDefX = 10;
715 fCanvasDefY = 10;
719 fPadBottomMargin= 0.1;
720 fPadTopMargin = 0.1;
721 fPadLeftMargin = 0.1;
722 fPadRightMargin = 0.1;
725 fPadTickX = 0;
726 fPadTickY = 0;
727 fFuncColor = 1;
728 fFuncStyle = 1;
729 fFuncWidth = 3;
730 fGridColor = 0;
731 fGridStyle = 3;
732 fGridWidth = 1;
733 fHatchesSpacing = 1;
735 fHistLineColor = 1;
736 fHistFillColor = 0;
737 fHistFillStyle = 1001;
738 fHistLineStyle = 1;
739 fHistLineWidth = 1;
741 fHistTopMargin = 0.05;
742 fFrameLineColor = 1;
743 fFrameFillColor = 0;
744 fFrameFillStyle = 1001;
745 fFrameLineStyle = 1;
746 fFrameLineWidth = 1;
749 fBarWidth = 1;
750 fBarOffset = 0;
753 fDrawBorder = 0;
754 fOptLogx = 0;
755 fOptLogy = 0;
756 fOptLogz = 0;
757 fOptDate = 0;
758 fOptFile = 0;
759 fOptFit = 0;
760 fOptStat = 1;
761 fOptTitle = 1;
762 fEndErrorSize = 2;
763 fErrorX = 0.5;
764 fScreenFactor = 1;
766 fStatTextColor = 1;
767 fStatBorderSize = 2;
768 fStatFont = 62;
769 fStatFontSize = 0;
770 fStatStyle = 1001;
771 fStatW = 0.20;
772 fStatH = 0.16;
773 fStatX = 0.98;
774 fStatY = 0.995;
776 SetFitFormat();
778 fTitleAlign = 13;
780 fTitleTextColor = 1;
781 fTitleFont = 62;
782 fTitleFontSize = 0;
783 fTitleStyle = 1001;
785 fTitleW = 0;
786 fTitleH = 0;
787 fTitleX = 0.01;
788 fTitleY = 0.995;
790 fShowEditor = 0;
791 fShowToolBar = 0;
792 fLegoInnerR = 0.5;
793 fHeaderPS = "";
794 fTitlePS = "";
796 fNumberContours = 20;
798 fLegendFont = 62;
799 fLegendTextSize = 0.,
801 fLegendFillStyle = 1001;
802 fImageScaling = 1.;
803
804 SetDateX();
805 SetDateY();
806 fAttDate.SetTextSize(0.025);
810 SetCapLinePS();
812 SetLineStyleString(1," ");
813 SetLineStyleString(2,"12 12");
814 SetLineStyleString(3,"4 8");
815 SetLineStyleString(4,"12 16 4 16");
816 SetLineStyleString(5,"20 12 4 12");
817 SetLineStyleString(6,"20 12 4 12 4 12 4 12");
818 SetLineStyleString(7,"20 20");
819 SetLineStyleString(8,"20 12 4 12 4 12");
820 SetLineStyleString(9,"80 20");
821 SetLineStyleString(10,"80 40 4 40");
822 for (Int_t i=11;i<30;i++) SetLineStyleString(i," ");
823
824 SetPaperSize();
825
826 SetPalette();
827
828 fTimeOffset = 788918400; // UTC time at 01/01/95
829
831 fCandleBoxRange = 0.5;
834
836
841 fAxisMaxDigits = 5;
842
843 TString style_name = opt;
844
845 if (strcmp(style_name,"Modern") == 0) {
846 // Modern style
852 SetPadColor(0);
853 SetStatColor(0);
854 SetTitleFont(42,"");
855 SetLabelFont(42,"x");
856 SetTitleFont(42,"x");
857 SetLabelFont(42,"y");
858 SetTitleFont(42,"y");
859 SetLabelFont(42,"z");
860 SetTitleFont(42,"z");
861 SetStatFont(42);
862 SetLabelSize(0.035,"x");
863 SetTitleSize(0.035,"x");
864 SetLabelSize(0.035,"y");
865 SetTitleSize(0.035,"y");
866 SetLabelSize(0.035,"z");
867 SetTitleSize(0.035,"z");
868 SetTitleSize(0.050,"");
869 SetTitleAlign(23);
870 SetTitleX(0.5);
873 SetTitleStyle(0);
874 SetTitleOffset(0.,"Y");
876 SetOptStat(1111);
877 SetStatY(0.935);
881 SetLegendFillStyle(1001);
882 SetLegendFont(42);
884 SetFuncWidth(2);
885 SetFuncColor(2);
886 }
887 if (strcmp(style_name,"Plain") == 0) {
891 SetPadColor(0);
895 SetStatColor(0);
898 return;
899 }
900 if (strcmp(style_name,"Bold") == 0) {
901 SetPalette(1,nullptr);
902 SetCanvasColor(10);
906 SetPadColor(10);
907 SetPadTickX(1);
908 SetPadTickY(1);
909 SetPadBottomMargin(0.15);
910 SetPadLeftMargin(0.15);
913 SetFuncWidth(3);
915 SetLineWidth(3);
916 SetLabelSize(0.05,"xyz");
917 SetLabelOffset(0.01,"y");
918 SetLabelColor(kBlue,"xy");
919 SetTitleSize(0.06,"xyz");
920 SetTitleOffset(1.3,"Y");
923 SetStatColor(10);
924 return;
925 }
926 if (strcmp(style_name,"Video") == 0) {
927 SetPalette(1,nullptr);
928 SetCanvasColor(10);
932 SetPadColor(10);
933 SetPadTickX(1);
934 SetPadTickY(1);
936 SetPadLeftMargin(0.2);
939 SetLabelSize(0.06,"xyz");
940 SetLabelColor(kBlue,"xyz");
941 SetTitleSize(0.08,"xyz");
944 SetStatColor(10);
945 SetFuncWidth(8);
947 SetLineWidth(3);
948 return;
949 }
950 if (strcmp(style_name,"Pub") == 0) {
951 SetOptTitle(0);
952 SetOptStat(0);
953 SetPalette(8,nullptr);
954 SetCanvasColor(10);
958 SetPadColor(10);
959 SetPadTickX(1);
960 SetPadTickY(1);
961 SetPadBottomMargin(0.15);
962 SetPadLeftMargin(0.15);
965 SetFuncWidth(3);
967 SetLineWidth(3);
968 SetLabelSize(0.05,"xyz");
969 SetLabelOffset(0.01,"y");
970 SetLabelColor(kBlack,"xyz");
971 SetTitleSize(0.06,"xyz");
972 SetTitleOffset(1.3,"y");
975 return;
976 }
977 if (strcmp(style_name,"ATLAS") == 0) {
983 SetPadColor(0);
984 SetStatColor(0);
985 SetPaperSize(20,26);
986 SetPadTopMargin(0.05);
987 SetPadRightMargin(0.05);
988 SetPadBottomMargin(0.16);
989 SetPadLeftMargin(0.16);
990 SetTitleXOffset(1.4);
991 SetTitleYOffset(1.4);
992 Int_t font = 42;
993 Double_t tsize=0.05;
994 SetTextFont(font);
996 SetLabelFont(font,"x");
997 SetTitleFont(font,"x");
998 SetLabelFont(font,"y");
999 SetTitleFont(font,"y");
1000 SetLabelFont(font,"z");
1001 SetTitleFont(font,"z");
1002 SetLabelSize(tsize,"x");
1003 SetTitleSize(tsize,"x");
1004 SetLabelSize(tsize,"y");
1005 SetTitleSize(tsize,"y");
1006 SetLabelSize(tsize,"z");
1007 SetTitleSize(tsize,"z");
1008 SetMarkerStyle(20);
1009 SetMarkerSize(1.2);
1010 SetHistLineWidth(2.);
1011 SetLineStyleString(2,"[12 12]");
1012 SetEndErrorSize(0.);
1013 SetOptTitle(0);
1014 SetOptStat(0);
1015 SetOptFit(0);
1016 SetPadTickX(1);
1017 SetPadTickY(1);
1018 return;
1019 }
1020 if (strcmp(style_name,"BELLE2") == 0) {
1021 Int_t icol=0;
1029 SetPaperSize(20,26);
1030 SetPadTopMargin(0.05);
1031 SetPadRightMargin(0.05);
1032 SetPadBottomMargin(0.16);
1033 SetPadLeftMargin(0.16);
1034 SetTitleXOffset(1.0);
1035 SetTitleYOffset(1.0);
1036 Int_t font=42;
1037 Double_t tsize=0.05;
1038 SetTextFont(font);
1040 SetLabelFont(font,"x");
1041 SetTitleFont(font,"x");
1042 SetLabelFont(font,"y");
1043 SetTitleFont(font,"y");
1044 SetLabelFont(font,"z");
1045 SetTitleFont(font,"z");
1046 SetLabelSize(tsize,"x");
1047 SetTitleSize(.065,"x");
1048 SetLabelSize(tsize,"y");
1049 SetTitleSize(.065,"y");
1050 SetLabelSize(tsize,"z");
1051 SetTitleSize(.065,"z");
1052 SetTitleOffset(1.1,"x");
1053 SetTitleOffset(1.1,"y");
1054 SetTitleOffset(1.1,"z");
1055 SetLabelOffset(0.015,"x");
1056 SetLabelOffset(0.015,"y");
1057 SetLabelOffset(0.015,"z");
1058 SetTickLength(0.03,"x");
1059 SetTickLength(0.02,"y");
1060 SetMarkerStyle(20);
1061 SetMarkerSize(0.9);
1063 SetLineStyleString(2,"[12 12]");
1064 SetErrorX(0.001);
1065 SetEndErrorSize(0.);
1066 SetOptTitle(0);
1067 SetOptStat(0);
1068 SetOptFit(0);
1069 SetPadTickX(0);
1070 SetPadTickY(0);
1072 }
1073}
1074
1075////////////////////////////////////////////////////////////////////////////////
1076/// Return number of divisions.
1077
1079{
1080 Int_t ax = AxisChoice(axis);
1081 if (ax == 1) return fXaxis.GetNdivisions();
1082 if (ax == 2) return fYaxis.GetNdivisions();
1083 if (ax == 3) return fZaxis.GetNdivisions();
1084 return 0;
1085}
1086
1087////////////////////////////////////////////////////////////////////////////////
1088/// Return the axis color number in the axis.
1089
1091{
1092 Int_t ax = AxisChoice(axis);
1093 if (ax == 1) return fXaxis.GetAxisColor();
1094 if (ax == 2) return fYaxis.GetAxisColor();
1095 if (ax == 3) return fZaxis.GetAxisColor();
1096 return 0;
1097}
1098
1099////////////////////////////////////////////////////////////////////////////////
1100/// Return color number i in current palette.
1101
1103{
1104 return TColor::GetColorPalette(i);
1105}
1106
1107////////////////////////////////////////////////////////////////////////////////
1108/// Return the label color number in the axis.
1109
1111{
1112 Int_t ax = AxisChoice(axis);
1113 if (ax == 1) return fXaxis.GetLabelColor();
1114 if (ax == 2) return fYaxis.GetLabelColor();
1115 if (ax == 3) return fZaxis.GetLabelColor();
1116 return 0;
1117}
1118
1119////////////////////////////////////////////////////////////////////////////////
1120/// Return label font.
1121
1123{
1124 Int_t ax = AxisChoice(axis);
1125 if (ax == 1) return fXaxis.GetLabelFont();
1126 if (ax == 2) return fYaxis.GetLabelFont();
1127 if (ax == 3) return fZaxis.GetLabelFont();
1128 return 0;
1129}
1130
1131////////////////////////////////////////////////////////////////////////////////
1132/// Return label offset.
1133
1135{
1136 Int_t ax = AxisChoice(axis);
1137 if (ax == 1) return fXaxis.GetLabelOffset();
1138 if (ax == 2) return fYaxis.GetLabelOffset();
1139 if (ax == 3) return fZaxis.GetLabelOffset();
1140 return 0;
1141}
1142
1143////////////////////////////////////////////////////////////////////////////////
1144/// Return label size.
1145
1147{
1148 Int_t ax = AxisChoice(axis);
1149 if (ax == 1) return fXaxis.GetLabelSize();
1150 if (ax == 2) return fYaxis.GetLabelSize();
1151 if (ax == 3) return fZaxis.GetLabelSize();
1152 return 0;
1153}
1154
1155////////////////////////////////////////////////////////////////////////////////
1156/// Method returns maximum number of digits permitted for the axis labels above which the
1157/// notation with 10^N is used. See @ref SetAxisMaxDigits for more details
1159{
1160 return fAxisMaxDigits;
1161}
1162
1163////////////////////////////////////////////////////////////////////////////////
1164/// Return line style string (used by PostScript).
1165/// See SetLineStyleString for more explanations
1166
1168{
1169 if (i < 1 || i > 29) return fLineStyle[0].Data();
1170 return fLineStyle[i].Data();
1171}
1172
1173////////////////////////////////////////////////////////////////////////////////
1174/// Return number of colors in the color palette.
1175
1180
1181////////////////////////////////////////////////////////////////////////////////
1182/// Set paper size for PostScript output.
1183
1189
1190////////////////////////////////////////////////////////////////////////////////
1191/// Return tick length.
1192
1194{
1195 Int_t ax = AxisChoice(axis);
1196 if (ax == 1) return fXaxis.GetTickLength();
1197 if (ax == 2) return fYaxis.GetTickLength();
1198 if (ax == 3) return fZaxis.GetTickLength();
1199 return 0;
1200}
1201
1202////////////////////////////////////////////////////////////////////////////////
1203/// Return title color.
1204
1206{
1207 Int_t ax = AxisChoice(axis);
1208 if (ax == 1) return fXaxis.GetTitleColor();
1209 if (ax == 2) return fYaxis.GetTitleColor();
1210 if (ax == 3) return fZaxis.GetTitleColor();
1211 return fTitleTextColor;
1212}
1213
1214////////////////////////////////////////////////////////////////////////////////
1215/// Return title font.
1216
1218{
1219 Int_t ax = AxisChoice(axis);
1220 if (ax == 1) return fXaxis.GetTitleFont();
1221 if (ax == 2) return fYaxis.GetTitleFont();
1222 if (ax == 3) return fZaxis.GetTitleFont();
1223 return fTitleFont;
1224}
1225
1226////////////////////////////////////////////////////////////////////////////////
1227/// Return title offset.
1228
1230{
1231 Int_t ax = AxisChoice(axis);
1232 if (ax == 1) return fXaxis.GetTitleOffset();
1233 if (ax == 2) return fYaxis.GetTitleOffset();
1234 if (ax == 3) return fZaxis.GetTitleOffset();
1235 return 0;
1236}
1237
1238////////////////////////////////////////////////////////////////////////////////
1239/// Return title size.
1240
1242{
1243 Int_t ax = AxisChoice(axis);
1244 if (ax == 1) return fXaxis.GetTitleSize();
1245 if (ax == 2) return fYaxis.GetTitleSize();
1246 if (ax == 3) return fZaxis.GetTitleSize();
1247 return fTitleFontSize;
1248}
1249
1250////////////////////////////////////////////////////////////////////////////////
1251/// Copy this style to gStyle.
1252
1254{
1255 Copy(*gStyle);
1256}
1257
1258////////////////////////////////////////////////////////////////////////////////
1259/// Define the color model used by TPostScript and TPDF (RGB or CMYK).
1260/// CMY and CMYK models are subtractive color models unlike RGB which is
1261/// additive. They are mainly used for printing purposes. CMY means Cyan Magenta
1262/// Yellow. To convert RGB to CMY it is enough to do: C=1-R, M=1-G and Y=1-B.
1263/// CMYK has one more component K (black). The conversion from RGB to CMYK is:
1264/// ~~~ {.cpp}
1265/// Double_t Black = TMath::Min(TMath::Min(1-Red,1-Green),1-Blue);
1266/// Double_t Cyan = (1-Red-Black)/(1-Black);
1267/// Double_t Magenta = (1-Green-Black)/(1-Black);
1268/// Double_t Yellow = (1-Blue-Black)/(1-Black);
1269/// ~~~
1270/// CMYK adds the black component which allows better quality for black
1271/// printing. PostScript and PDF support the CMYK model.
1272///
1273/// - c = 0 means TPostScript and TPDF will use RGB color model (default)
1274/// - c = 1 means TPostScript and TPDF will use CMYK color model
1275
1280
1281////////////////////////////////////////////////////////////////////////////////
1282/// If the argument zero=kTRUE the minimum value for the Y axis of 1-d histograms
1283/// is set to 0.
1284///
1285/// If the minimum bin content is greater than 0 and TH1::SetMinimum
1286/// has not been called.
1287/// Otherwise the minimum is based on the minimum bin content.
1288
1293
1294////////////////////////////////////////////////////////////////////////////////
1295/// Set the number of divisions to draw an axis.
1296/// ndiv : Number of divisions.
1297/// ~~~ {.cpp}
1298/// n = N1 + 100*N2 + 10000*N3
1299/// N1=number of primary divisions.
1300/// N2=number of secondary divisions.
1301/// N3=number of 3rd divisions.
1302/// e.g.:
1303/// nndi=0 --> no tick marks.
1304/// nndi=2 --> 2 divisions, one tick mark in the middle
1305/// of the axis.
1306/// ~~~
1307/// axis specifies which axis ("x","y","z"), default = "x"
1308/// if axis="xyz" set all 3 axes
1309
1311{
1312 TString opt = axis;
1313 opt.ToLower();
1314 if (opt.Contains("x")) fXaxis.SetNdivisions(n);
1315 if (opt.Contains("y")) fYaxis.SetNdivisions(n);
1316 if (opt.Contains("z")) fZaxis.SetNdivisions(n);
1317}
1318
1319////////////////////////////////////////////////////////////////////////////////
1320/// Set color to draw the axis line and tick marks.
1321/// axis specifies which axis ("x","y","z"), default = "x"
1322/// if axis="xyz" set all 3 axes
1323
1325{
1326 TString opt = axis;
1327 opt.ToLower();
1328
1329 if (opt.Contains("x")) fXaxis.SetAxisColor(color);
1330 if (opt.Contains("y")) fYaxis.SetAxisColor(color);
1331 if (opt.Contains("z")) fZaxis.SetAxisColor(color);
1332}
1333
1334////////////////////////////////////////////////////////////////////////////////
1335/// Set the size (in pixels) of the small lines drawn at the
1336/// end of the error bars (TH1 or TGraphErrors).
1337///
1338/// The default value is 2 pixels.
1339/// Set np=0 to remove these lines
1340
1342{
1343 if (np >= 0) fEndErrorSize = np;
1344 else fEndErrorSize = 0;
1345}
1346
1347////////////////////////////////////////////////////////////////////////////////
1348/// Define a string to be inserted in the Postscript header.
1349///
1350/// The string in header will be added to the Postscript file
1351/// immediately following the %%Page line
1352/// For example, this string may contain special Postscript instructions like
1353/// ~~~ {.cpp}
1354/// 200 200 translate
1355/// ~~~
1356/// the following header string will print the string "my annotation" at the
1357/// bottom left corner of the page (outside the user area)
1358/// ~~~ {.cpp}
1359/// "gsave 100 -100 t 0 r 0 0 m /Helvetica-Bold findfont 56 sf 0 0 m ( my annotation ) show gr"
1360/// ~~~
1361/// This information is used in TPostScript::Initialize
1362
1363void TStyle::SetHeaderPS(const char *header)
1364{
1365 fHeaderPS = header;
1366}
1367
1368////////////////////////////////////////////////////////////////////////////////
1369/// Sets the `fIsReading` member to reading (default=kTRUE).
1370///
1371/// `fIsReading` (used via `gStyle->IsReading()`) can be used in
1372/// the functions `myclass::UseCurrentStyle` to read from the current style
1373/// or write to the current style
1374
1379
1380////////////////////////////////////////////////////////////////////////////////
1381/// Define a string to be used in the %%Title of the Postscript files.
1382/// If this string is not defined, ROOT will use the canvas title.
1383
1385{
1386 fTitlePS = pstitle;
1387}
1388
1389////////////////////////////////////////////////////////////////////////////////
1390/// Set axis labels color.
1391/// axis specifies which axis ("x","y","z"), default = "x"
1392/// if axis="xyz" set all 3 axes
1393
1395{
1396 TString opt = axis;
1397 opt.ToLower();
1398
1399 if (opt.Contains("x")) fXaxis.SetLabelColor(color);
1400 if (opt.Contains("y")) fYaxis.SetLabelColor(color);
1401 if (opt.Contains("z")) fZaxis.SetLabelColor(color);
1402}
1403
1404////////////////////////////////////////////////////////////////////////////////
1405/// Set font number used to draw axis labels.
1406/// - font : Text font code = 10*fontnumber + precision
1407/// - Font numbers must be between 1 and 14
1408/// - precision = 1 fast hardware fonts (steps in the size)
1409/// - precision = 2 scalable and rotatable hardware fonts
1410/// The default font number is 62.
1411/// axis specifies which axis ("x","y","z"), default = "x"
1412/// if axis="xyz" set all 3 axes
1413
1415{
1416 TString opt = axis;
1417 opt.ToLower();
1418
1419 if (opt.Contains("x")) fXaxis.SetLabelFont(font);
1420 if (opt.Contains("y")) fYaxis.SetLabelFont(font);
1421 if (opt.Contains("z")) fZaxis.SetLabelFont(font);
1422}
1423
1424////////////////////////////////////////////////////////////////////////////////
1425/// Set offset between axis and axis labels.
1426/// The offset is expressed as a percent of the pad height.
1427/// axis specifies which axis ("x","y","z"), default = "x"
1428/// if axis="xyz" set all 3 axes
1429
1431{
1432 TString opt = axis;
1433 opt.ToLower();
1434
1435 if (opt.Contains("x")) fXaxis.SetLabelOffset(offset);
1436 if (opt.Contains("y")) fYaxis.SetLabelOffset(offset);
1437 if (opt.Contains("z")) fZaxis.SetLabelOffset(offset);
1438}
1439
1440////////////////////////////////////////////////////////////////////////////////
1441/// Set size of axis labels. The size is expressed as a percent of the pad height.
1442/// axis specifies which axis ("x","y","z"), default = "x"
1443/// if axis="xyz" set all 3 axes
1444
1446{
1447 TString opt = axis;
1448 opt.ToLower();
1449
1450 if (opt.Contains("x")) fXaxis.SetLabelSize(size);
1451 if (opt.Contains("y")) fYaxis.SetLabelSize(size);
1452 if (opt.Contains("z")) fZaxis.SetLabelSize(size);
1453}
1454
1455////////////////////////////////////////////////////////////////////////////////
1456/// Set line style string using the PostScript convention.
1457/// A line is a suite of segments, each segment is described by the number of
1458/// pixels. The initial and alternating elements (second, fourth, and so on)
1459/// are the dashes, and the others spaces between dashes.
1460///
1461/// Default fixed line styles are pre-defined as:
1462/// ~~~ {.cpp}
1463/// linestyle 1 "[]" solid
1464/// linestyle 2 "[12 12]" dashed
1465/// linestyle 3 "[4 8]" dotted
1466/// linestyle 4 "[12 16 4 16]" dash-dotted
1467/// ~~~
1468/// For example the following lines define the line style 5 to 9.
1469/// ~~~ {.cpp}
1470/// gStyle->SetLineStyleString(5,"20 12 4 12");
1471/// gStyle->SetLineStyleString(6,"20 12 4 12 4 12 4 12");
1472/// gStyle->SetLineStyleString(7,"20 20");
1473/// gStyle->SetLineStyleString(8,"20 12 4 12 4 12");
1474/// gStyle->SetLineStyleString(9,"80 20");
1475/// ~~~
1476/// \image html base_linestyle.png
1477/// Note:
1478/// - Up to 30 different styles may be defined.
1479/// - The opening and closing brackets may be omitted
1480/// - It is recommended to use 4 as the smallest segment length and multiple of
1481/// 4 for other lengths.
1482/// - The line style 1 to 10 are predefined. 1 to 4 cannot be changed.
1483
1485{
1486 if (!text) text = "";
1487 char *l;
1488 Int_t nch = strlen(text);
1489 char *st = new char[nch+10];
1490 snprintf(st,nch+10," ");
1491 strlcat(st,text,nch+10);
1492 l = strstr(st,"["); if (l) l[0] = ' ';
1493 l = strstr(st,"]"); if (l) l[0] = ' ';
1494 if (i >= 1 && i <= 29) fLineStyle[i] = st;
1495 delete [] st;
1496}
1497
1498////////////////////////////////////////////////////////////////////////////////
1499/// Set the default number of contour levels when drawing 2-d plots.
1500
1502{
1503 if (number > 0 && number < 1000) {
1504 fNumberContours = number;
1505 return;
1506 }
1507
1508 Error("SetNumberContours","Illegal number of contours: %d, must be > 0 and < 1000",number);
1509}
1510
1511////////////////////////////////////////////////////////////////////////////////
1512/// If optdate is non null, the current date/time will be printed in the canvas.
1513/// The position of the date string can be controlled by:
1514/// optdate = 10*format + mode
1515/// - mode = 1 (default) date is printed in the bottom/left corner.
1516/// - mode = 2 date is printed in the bottom/right corner.
1517/// - mode = 3 date is printed in the top/right corner.
1518/// - format = 0 (default) date has the format like: "Wed Sep 25 17:10:35 2002"
1519/// - format = 1 date has the format like: "2002-09-25"
1520/// - format = 2 date has the format like: "2002-09-25 17:10:35"
1521///
1522/// examples:
1523/// - optdate = 1 date like "Wed Sep 25 17:10:35 2002" in the bottom/left corner.
1524/// - optdate = 13 date like "2002-09-25" in the top/right corner.
1525///
1526/// The date position can also be controlled by:
1527/// gStyle->SetDateX(x); x in NDC
1528/// gStyle->SetDateY(y); y in NDC
1529///
1530/// The date text attributes can be changed with:
1531/// ~~~ {.cpp}
1532/// gStyle->GetAttDate()->SetTextFont(font=62);
1533/// gStyle->GetAttDate()->SetTextSize(size=0.025);
1534/// gStyle->GetAttDate()->SetTextAngle(angle=0);
1535/// gStyle->GetAttDate()->SetTextAlign(align=11);
1536/// gStyle->GetAttDate()->SetTextColor(color=1);
1537/// ~~~
1538/// The current date attributes can be obtained via:
1539/// ~~~ {.cpp}
1540/// gStyle->GetAttDate()->GetTextxxxx();
1541/// ~~~
1542/// When the date option is active, a text object is created when the pad
1543/// paint its list of primitives. The text object is named "DATE".
1544/// The DATE attributes can also be edited interactively (position
1545/// and attributes) via the normal context menu.
1546
1548{
1549 fOptDate = optdate;
1550 Int_t mode = optdate % 10;
1551 if (mode == 1) {
1552 SetDateX(0.01);
1553 SetDateY(0.01);
1555 }
1556 if (mode == 2) {
1557 SetDateX(0.99);
1558 SetDateY(0.01);
1560 }
1561 if (mode == 3) {
1562 SetDateX(0.99);
1563 SetDateY(0.99);
1565 }
1566}
1567
1568////////////////////////////////////////////////////////////////////////////////
1569/// The type of information about fit parameters printed in the histogram
1570/// statistics box can be selected via the parameter `mode`.
1571/// The parameter mode can be = `pcev`:
1572/// - p = 1; print Probability
1573/// - c = 1; print Chisquare/Number of degrees of freedom
1574/// - e = 1; print errors (if e=1, v must be 1)
1575/// - v = 1; print name/values of parameters
1576/// Example: `gStyle->SetOptFit(1011);`
1577/// print fit probability, parameter names/values and errors.
1578/// - When "v"=1 is specified, only the non-fixed parameters are shown.
1579/// - When "v"=2 all parameters are shown.
1580///
1581/// #### Notes:
1582///
1583/// - never call `SetOptFit(000111);` but `SetOptFit(111)`, 000111 will
1584/// be taken as an octal number !!
1585/// - `gStyle->SetOptFit(1)` is a shortcut allowing to set the most common
1586/// case and is equivalent to `gStyle->SetOptFit(111)`
1587/// - At ROOT startup the option fit is set to `0`. So, to see the fit parameters
1588/// on all plot resulting from a fit, a call to `gStyle->SetOptFit()` with a
1589/// non null value should be done. One can put it in the `rootlogon.C` file to
1590/// always have it.
1591///
1592/// see also SetOptStat below.
1593
1595{
1596 fOptFit = mode;
1597 if (gPad) {
1598 TIter next(gPad->GetListOfPrimitives());
1599 while (auto obj = next()) {
1600 TObject *stats = obj->FindObject("stats");
1601 if (stats) stats->SetBit(kTakeStyle);
1602 }
1603 gPad->Modified(); gPad->Update();
1604 }
1605}
1606
1607////////////////////////////////////////////////////////////////////////////////
1608/// The type of information printed in the histogram statistics box
1609/// can be selected via the parameter mode.
1610/// The parameter mode can be = `ksiourmen`
1611/// - k = 1; kurtosis printed
1612/// - k = 2; kurtosis and kurtosis error printed
1613/// - s = 1; skewness printed
1614/// - s = 2; skewness and skewness error printed
1615/// - i = 1; integral of bins printed
1616/// - i = 2; integral of bins with option "width" printed
1617/// - o = 1; number of overflows printed
1618/// - u = 1; number of underflows printed
1619/// - r = 1; rms printed
1620/// - r = 2; rms and rms error printed
1621/// - m = 1; mean value printed
1622/// - m = 2; mean and mean error values printed
1623/// - e = 1; number of entries printed
1624/// - n = 1; name of histogram is printed
1625///
1626/// Example: `gStyle->SetOptStat(11);`
1627/// print only name of histogram and number of entries.
1628/// `gStyle->SetOptStat(1101);` displays the name of histogram, mean value and RMS.
1629///
1630/// #### Notes:
1631///
1632/// - never call `SetOptStat(000111);` but `SetOptStat(111)`, 000111 will
1633/// be taken as an octal number !!
1634/// - `SetOptStat(1)` is s shortcut allowing to set the most common case, and is
1635/// taken as `SetOptStat(1111)` (for backward compatibility with older versions.
1636/// If you want to print only the name of the histogram call `SetOptStat(1000000001)`.
1637/// - that in case of 2-D histograms, when selecting just underflow (10000)
1638/// or overflow (100000), the stats box will show all combinations
1639/// of underflow/overflows and not just one single number!
1640
1642{
1643 fOptStat = mode;
1644 if (gPad) {
1645 TIter next(gPad->GetListOfPrimitives());
1646 while (auto obj = next()) {
1647 TObject *stats = obj->FindObject("stats");
1648 if (stats) stats->SetBit(kTakeStyle);
1649 }
1650 gPad->Modified(); gPad->Update();
1651 }
1652}
1653
1654////////////////////////////////////////////////////////////////////////////////
1655/// The parameter mode can be any combination of kKsSiourRmMen
1656/// - k : kurtosis printed
1657/// - K : kurtosis and kurtosis error printed
1658/// - s : skewness printed
1659/// - S : skewness and skewness error printed
1660/// - i : integral of bins printed
1661/// - I : integral of bins with option "width" printed
1662/// - o : number of overflows printed
1663/// - u : number of underflows printed
1664/// - r : rms printed
1665/// - R : rms and rms error printed
1666/// - m : mean value printed
1667/// - M : mean value mean error values printed
1668/// - e : number of entries printed
1669/// - n : name of histogram is printed
1670///
1671/// Example: `gStyle->SetOptStat("ne");`
1672/// print only name of histogram and number of entries.
1673///
1674/// - `gStyle->SetOptStat("n")` print only the name of the histogram
1675/// - `gStyle->SetOptStat("nemr")` is the default
1676
1678{
1679 Int_t mode = 0;
1680
1681 TString opt = stat;
1682
1683 if (opt.Contains("n")) mode+=1;
1684 if (opt.Contains("e")) mode+=10;
1685 if (opt.Contains("m")) mode+=100;
1686 if (opt.Contains("M")) mode+=200;
1687 if (opt.Contains("r")) mode+=1000;
1688 if (opt.Contains("R")) mode+=2000;
1689 if (opt.Contains("u")) mode+=10000;
1690 if (opt.Contains("o")) mode+=100000;
1691 if (opt.Contains("i")) mode+=1000000;
1692 if (opt.Contains("I")) mode+=2000000;
1693 if (opt.Contains("s")) mode+=10000000;
1694 if (opt.Contains("S")) mode+=20000000;
1695 if (opt.Contains("k")) mode+=100000000;
1696 if (opt.Contains("K")) mode+=200000000;
1697 if (mode == 1) mode = 1000000001;
1698
1700}
1701
1702////////////////////////////////////////////////////////////////////////////////
1703/// Set paper size for PostScript output.
1704
1706{
1707 switch (size) {
1708 case kA4:
1709 SetPaperSize(20, 26);
1710 break;
1711 case kUSLetter:
1712 SetPaperSize(20, 24);
1713 break;
1714 default:
1715 Error("SetPaperSize", "illegal paper size %d\n", (int)size);
1716 break;
1717 }
1718}
1719
1720////////////////////////////////////////////////////////////////////////////////
1721/// Set paper size for PostScript output.
1722/// The paper size is specified in centimeters. Default is 20x26.
1723/// See also TPad::Print
1724
1730
1731////////////////////////////////////////////////////////////////////////////////
1732/// Set the tick marks length for an axis.
1733/// axis specifies which axis ("x","y","z"), default = "x"
1734/// if axis="xyz" set all 3 axes
1735
1737{
1738 TString opt = axis;
1739 opt.ToLower();
1740
1741 if (opt.Contains("x")) fXaxis.SetTickLength(length);
1742 if (opt.Contains("y")) fYaxis.SetTickLength(length);
1743 if (opt.Contains("z")) fZaxis.SetTickLength(length);
1744}
1745
1746////////////////////////////////////////////////////////////////////////////////
1747/// - if axis =="x" set the X axis title color
1748/// - if axis =="y" set the Y axis title color
1749/// - if axis =="z" set the Z axis title color
1750///
1751/// any other value of axis will set the pad title color
1752///
1753/// if axis="xyz" set all 3 axes
1754
1756{
1757 TString opt = axis;
1758 opt.ToLower();
1759
1760 Bool_t set = kFALSE;
1761 if (opt.Contains("x")) { fXaxis.SetTitleColor(color); set = kTRUE; }
1762 if (opt.Contains("y")) { fYaxis.SetTitleColor(color); set = kTRUE; }
1763 if (opt.Contains("z")) { fZaxis.SetTitleColor(color); set = kTRUE; }
1764 if (!set) fTitleColor = color;
1765}
1766
1767////////////////////////////////////////////////////////////////////////////////
1768/// - if axis =="x" set the X axis title font
1769/// - if axis =="y" set the Y axis title font
1770/// - if axis =="z" set the Z axis title font
1771///
1772/// any other value of axis will set the pad title font
1773///
1774/// if axis="xyz" set all 3 axes
1775
1777{
1778 TString opt = axis;
1779 opt.ToLower();
1780
1781 Bool_t set = kFALSE;
1782 if (opt.Contains("x")) { fXaxis.SetTitleFont(font); set = kTRUE; }
1783 if (opt.Contains("y")) { fYaxis.SetTitleFont(font); set = kTRUE; }
1784 if (opt.Contains("z")) { fZaxis.SetTitleFont(font); set = kTRUE; }
1785 if (!set) fTitleFont = font;
1786}
1787
1788////////////////////////////////////////////////////////////////////////////////
1789/// Specify a parameter offset to control the distance between the axis
1790/// and the axis title.
1791///
1792/// - offset = 1 means : use the default distance
1793/// - offset = 1.2 means: the distance will be 1.2*(default distance)
1794/// - offset = 0.8 means: the distance will be 0.8*(default distance)
1795///
1796/// axis specifies which axis ("x","y","z"), default = "x"
1797/// if axis="xyz" set all 3 axes
1798
1800{
1801 TString opt = axis;
1802 opt.ToLower();
1803
1804 if (opt.Contains("x")) fXaxis.SetTitleOffset(offset);
1805 if (opt.Contains("y")) fYaxis.SetTitleOffset(offset);
1806 if (opt.Contains("z")) fZaxis.SetTitleOffset(offset);
1807}
1808
1809////////////////////////////////////////////////////////////////////////////////
1810/// - if axis =="x" set the X axis title size
1811/// - if axis =="y" set the Y axis title size
1812/// - if axis =="z" set the Z axis title size
1813///
1814/// any other value of axis will set the pad title size
1815///
1816/// if axis="xyz" set all 3 axes
1817
1819{
1820 TString opt = axis;
1821 opt.ToLower();
1822
1823 Bool_t set = kFALSE;
1824 if (opt.Contains("x")) { fXaxis.SetTitleSize(size); set = kTRUE; }
1825 if (opt.Contains("y")) { fYaxis.SetTitleSize(size); set = kTRUE; }
1826 if (opt.Contains("z")) { fZaxis.SetTitleSize(size); set = kTRUE; }
1827 if (!set) fTitleFontSize = size;
1828}
1829
1830////////////////////////////////////////////////////////////////////////////////
1831/// Method set X and Y offset of the axis 10^n notation.
1832/// It applies on axis belonging to an histogram (TAxis). It has no effect on standalone TGaxis.
1833/// It is in % of the pad size. It can be negative.
1834/// axis specifies which axis ("x","y"), default = "x"
1835/// if axis="xz" set the two axes
1836
1838{
1839 TString opt = axis;
1840 opt.ToLower();
1841
1842 if (opt.Contains("x")) {
1845 }
1846 if (opt.Contains("y")) {
1849 }
1850}
1851
1852////////////////////////////////////////////////////////////////////////////////
1853/// Method returns X and Y offset of the axis 10^n notation.
1854/// It applies on axis belonging to an histogram (TAxis)
1855
1857{
1858 TString opt = axis;
1859 opt.ToLower();
1860
1861 if (opt.Contains("x")) {
1864 } else if (opt.Contains("y")) {
1867 } else {
1868 xoff = yoff = 0.;
1869 }
1870}
1871
1872////////////////////////////////////////////////////////////////////////////////
1873/// Method set maximum number of digits permitted for the axis labels above which the
1874/// notation with 10^N is used. For example, to accept 6 digits number like 900000
1875/// on an axis call `gStyle->SetAxisMaxDigits(6)`. The default value is 5.
1876/// Warning: this function changes the max number of digits in all axes.
1877/// If you only want to change the digits of the current TGaxis instance, use
1878/// axis->SetNdivisions(N*1000000 + (axis->GetNdiv()%1000000))
1879/// instead of gStyle->SetAxisMaxDigits(N).
1880
1882{
1883 fAxisMaxDigits = (maxd > 1) ? maxd : 1;
1884}
1885
1886////////////////////////////////////////////////////////////////////////////////
1887/// See TColor::SetPalette.
1888
1890{
1891 TColor::SetPalette(ncolors,colors,alpha);
1892}
1893
1894////////////////////////////////////////////////////////////////////////////////
1895/// \see TColor::CreateColorTableFromFile, TColor::SetPalette
1897{
1898 TColor::CreateColorTableFromFile(fileName, alpha);
1899}
1900
1901////////////////////////////////////////////////////////////////////////////////
1902/// Change the time offset for time plotting.
1903/// Times are expressed in seconds. The corresponding numbers usually have 9
1904/// digits (or more if one takes into account fractions of seconds).
1905/// Thus, since it is very inconvenient to plot very large numbers on a scale,
1906/// one has to set an offset time that will be added to the axis beginning,
1907/// in order to plot times correctly and conveniently. A convenient way to
1908/// set the time offset is to use TDatime::Convert().
1909///
1910/// By default the time offset is set to 788918400 which corresponds to
1911/// 01/01/1995. This allows to have valid dates until 2072. The standard
1912/// UNIX time offset in 1970 allows only valid dates until 2030.
1913
1918
1919////////////////////////////////////////////////////////////////////////////////
1920/// Set option to strip decimals when drawing axis labels.
1921/// By default, TGaxis::PaintAxis removes trailing 0s after a dot
1922/// in the axis labels. Ex: {0,0.5,1,1.5,2,2.5, etc}
1923/// If this function is called with strip=kFALSE, TGAxis::PaintAxis will
1924/// draw labels with the same number of digits after the dot
1925/// Ex: (0.0,0.5,1.0,1.5,2.0,2.5,etc}
1926
1931
1932////////////////////////////////////////////////////////////////////////////////
1933/// By setting whisker-range for candle plot, one can force
1934/// the whiskers to cover the fraction of the distribution.
1935/// Set wRange between 0 and 1. Default is 1
1936/// gStyle->SetCandleWhiskerRange(0.95) will set all candle-charts to cover 95% of
1937/// the distribution with the whiskers.
1938/// Can only be used with the standard-whisker definition
1939
1941{
1942 if (wRange < 0)
1944 else if (wRange > 1)
1946 else
1948}
1949
1950////////////////////////////////////////////////////////////////////////////////
1951/// By setting box-range for candle plot, one can force the
1952/// box of the candle-chart to cover that given fraction of the distribution.
1953/// Set bRange between 0 and 1. Default is 0.5
1954/// gStyle->SetCandleBoxRange(0.68) will set all candle-charts to cover 68% of the
1955/// distribution by the box
1956
1958{
1959 if (bRange < 0)
1960 fCandleBoxRange = 0;
1961 else if (bRange > 1)
1962 fCandleBoxRange = 1;
1963 else
1965}
1966
1967////////////////////////////////////////////////////////////////////////////////
1968/// Set the line width of the circle marker of a candle plot ([1,5]).
1969
1971{
1973 Error("SetCandleCircleLineWidth", "illegal line width %d. It must be in the range [1,5]\n", (int)CircleLineWidth);
1975 return;
1976 }
1978}
1979
1980////////////////////////////////////////////////////////////////////////////////
1981/// Set the line width of the cross marker of a candle plot ([1,5]).
1982
1984{
1986 Error("SetCandleCrossLineWidth", "illegal line width %d. It must be in the range [1,5]\n", (int)CrossLineWidth);
1988 return;
1989 }
1991}
1992
1993
1994////////////////////////////////////////////////////////////////////////////////
1995/// Save the current style in a C++ macro file.
1996
1998{
1999 // Opens a file named filename or "Rootstyl.C"
2000 TString ff = filename && *filename ? filename : "Rootstyl.C";
2001
2002 // Computes the main method name.
2004 auto pos = fname.First('.');
2005
2006 if (pos == kNPOS) {
2007 sname = fname;
2008 ff.Append(".C");
2009 } else
2010 sname = fname(0, pos);
2011
2012 // Tries to open the file.
2013 std::ofstream out;
2014 out.open(ff.Data(), std::ios::out);
2015 if (!out.good()) {
2016 Error("SaveSource", "cannot open file: %s", ff.Data());
2017 return;
2018 }
2019
2020 // Writes macro header, date/time stamp as string, and the used Root version
2021 out << "// Macro generated from application: " << gApplication->Argv(0) << "\n";
2022 out << "// By ROOT version " << gROOT->GetVersion() << " on " << TDatime().AsSQLString() << "\n\n";
2023
2024 // Writes include.
2025 out << "#include \"TStyle.h\"\n\n";
2026
2027 // Writes the macro entry point equal to the fname
2028 out << "void " << sname << "()\n";
2029 out << "{\n";
2030
2032
2033 out << "}\n";
2034 out.close();
2035
2036 printf(" C++ macro file %s has been generated\n", gSystem->BaseName(ff));
2037}
2038
2039////////////////////////////////////////////////////////////////////////////////
2040/// Save primitive as a C++ statement(s) on output stream out
2041
2042void TStyle::SavePrimitive(std::ostream &out, Option_t * /*= ""*/)
2043{
2044 auto asBool = [](bool flag) { return flag ? "kTRUE" : "kFALSE"; };
2045
2046 out << " // Add the saved style to the current ROOT session.\n";
2047 out << " \n";
2048 out << " delete gROOT->GetStyle(\"" << GetName() << "\");\n";
2049
2051 out, Class(), "tmpStyle",
2052 TString::Format("\"%s\", \"%s\"", GetName(), TString(GetTitle()).ReplaceSpecialCppChars().Data()));
2053
2054 const char *prefix = " tmpStyle->";
2055
2056 // fXAxis, fYAxis and fZAxis
2057 out << prefix << "SetNdivisions(" << GetNdivisions("x") << ", \"x\");\n";
2058 out << prefix << "SetNdivisions(" << GetNdivisions("y") << ", \"y\");\n";
2059 out << prefix << "SetNdivisions(" << GetNdivisions("z") << ", \"z\");\n";
2060 out << prefix << "SetAxisColor(" << TColor::SavePrimitiveColor(GetAxisColor("x")) << ", \"x\");\n";
2061 out << prefix << "SetAxisColor(" << TColor::SavePrimitiveColor(GetAxisColor("y")) << ", \"y\");\n";
2062 out << prefix << "SetAxisColor(" << TColor::SavePrimitiveColor(GetAxisColor("z")) << ", \"z\");\n";
2063 out << prefix << "SetLabelColor(" << TColor::SavePrimitiveColor(GetLabelColor("x")) << ", \"x\");\n";
2064 out << prefix << "SetLabelColor(" << TColor::SavePrimitiveColor(GetLabelColor("y")) << ", \"y\");\n";
2065 out << prefix << "SetLabelColor(" << TColor::SavePrimitiveColor(GetLabelColor("z")) << ", \"z\");\n";
2066 out << prefix << "SetLabelFont(" << GetLabelFont("x") << ", \"x\");\n";
2067 out << prefix << "SetLabelFont(" << GetLabelFont("y") << ", \"y\");\n";
2068 out << prefix << "SetLabelFont(" << GetLabelFont("z") << ", \"z\");\n";
2069 out << prefix << "SetLabelOffset(" << GetLabelOffset("x") << ", \"x\");\n";
2070 out << prefix << "SetLabelOffset(" << GetLabelOffset("y") << ", \"y\");\n";
2071 out << prefix << "SetLabelOffset(" << GetLabelOffset("z") << ", \"z\");\n";
2072 out << prefix << "SetLabelSize(" << GetLabelSize("x") << ", \"x\");\n";
2073 out << prefix << "SetLabelSize(" << GetLabelSize("y") << ", \"y\");\n";
2074 out << prefix << "SetLabelSize(" << GetLabelSize("z") << ", \"z\");\n";
2075 out << prefix << "SetTickLength(" << GetTickLength("x") << ", \"x\");\n";
2076 out << prefix << "SetTickLength(" << GetTickLength("y") << ", \"y\");\n";
2077 out << prefix << "SetTickLength(" << GetTickLength("z") << ", \"z\");\n";
2078 out << prefix << "SetTitleOffset(" << GetTitleOffset("x") << ", \"x\");\n";
2079 out << prefix << "SetTitleOffset(" << GetTitleOffset("y") << ", \"y\");\n";
2080 out << prefix << "SetTitleOffset(" << GetTitleOffset("z") << ", \"z\");\n";
2081 out << prefix << "SetTitleSize(" << GetTitleSize("x") << ", \"x\");\n";
2082 out << prefix << "SetTitleSize(" << GetTitleSize("y") << ", \"y\");\n";
2083 out << prefix << "SetTitleSize(" << GetTitleSize("z") << ", \"z\");\n";
2084 out << prefix << "SetTitleColor(" << TColor::SavePrimitiveColor(GetTitleColor("x")) << ", \"x\");\n";
2085 out << prefix << "SetTitleColor(" << TColor::SavePrimitiveColor(GetTitleColor("y")) << ", \"y\");\n";
2086 out << prefix << "SetTitleColor(" << TColor::SavePrimitiveColor(GetTitleColor("z")) << ", \"z\");\n";
2087 out << prefix << "SetTitleFont(" << GetTitleFont("x") << ", \"x\");\n";
2088 out << prefix << "SetTitleFont(" << GetTitleFont("y") << ", \"y\");\n";
2089 out << prefix << "SetTitleFont(" << GetTitleFont("z") << ", \"z\");\n";
2090
2091 out << prefix << "SetExponentOffset(" << fXAxisExpXOffset << ", " << fXAxisExpYOffset << ", \"x\");\n";
2092 out << prefix << "SetExponentOffset(" << fYAxisExpXOffset << ", " << fYAxisExpYOffset << ", \"y\");\n";
2093 out << prefix << "SetAxisMaxDigits(" << GetAxisMaxDigits() << ");\n";
2094
2095 out << prefix << "SetBarWidth(" << GetBarWidth() << ");\n";
2096 out << prefix << "SetBarOffset(" << GetBarOffset() << ");\n";
2097 out << prefix << "SetDrawBorder(" << GetDrawBorder() << ");\n";
2098 out << prefix << "SetOptLogx(" << GetOptLogx() << ");\n";
2099 out << prefix << "SetOptLogy(" << GetOptLogy() << ");\n";
2100 out << prefix << "SetOptLogz(" << GetOptLogz() << ");\n";
2101 out << prefix << "SetOptDate(" << GetOptDate() << ");\n";
2102 out << prefix << "SetOptStat(" << GetOptStat() << ");\n";
2103 out << prefix << "SetOptTitle(" << GetOptTitle() << ");\n";
2104 out << prefix << "SetOptFit(" << GetOptFit() << ");\n";
2105 out << prefix << "SetNumberContours(" << GetNumberContours() << ");\n";
2106
2107 // fAttDate
2108 out << prefix << "GetAttDate()->SetTextFont(" << GetAttDate()->GetTextFont() << ");\n";
2109 out << prefix << "GetAttDate()->SetTextSize(" << GetAttDate()->GetTextSize() << ");\n";
2110 out << prefix << "GetAttDate()->SetTextAngle(" << GetAttDate()->GetTextAngle() << ");\n";
2111 out << prefix << "GetAttDate()->SetTextAlign(" << GetAttDate()->GetTextAlign() << ");\n";
2112 out << prefix << "GetAttDate()->SetTextColor(" << TColor::SavePrimitiveColor(GetAttDate()->GetTextColor()) << ");\n";
2113
2114 out << prefix << "SetDateX(" << GetDateX() << ");\n";
2115 out << prefix << "SetDateY(" << GetDateY() << ");\n";
2116 out << prefix << "SetEndErrorSize(" << GetEndErrorSize() << ");\n";
2117 out << prefix << "SetErrorX(" << GetErrorX() << ");\n";
2118 out << prefix << "SetFuncColor(" << TColor::SavePrimitiveColor(GetFuncColor()) << ");\n";
2119 out << prefix << "SetFuncStyle(" << GetFuncStyle() << ");\n";
2120 out << prefix << "SetFuncWidth(" << GetFuncWidth() << ");\n";
2121 out << prefix << "SetGridColor(" << TColor::SavePrimitiveColor(GetGridColor()) << ");\n";
2122 out << prefix << "SetGridStyle(" << GetGridStyle() << ");\n";
2123 out << prefix << "SetGridWidth(" << GetGridWidth() << ");\n";
2124 out << prefix << "SetLegendBorderSize(" << GetLegendBorderSize() << ");\n";
2125 out << prefix << "SetLegendFillColor(" << TColor::SavePrimitiveColor(GetLegendFillColor()) << ");\n";
2126 out << prefix << "SetLegendFillStyle(" << GetLegendFillStyle() << ");\n";
2127 out << prefix << "SetLegendFont(" << GetLegendFont() << ");\n";
2128 out << prefix << "SetLegendTextSize(" << GetLegendTextSize() << ");\n";
2129 out << prefix << "SetHatchesLineWidth(" << GetHatchesLineWidth() << ");\n";
2130 out << prefix << "SetHatchesSpacing(" << GetHatchesSpacing() << ");\n";
2131 out << prefix << "SetFrameFillColor(" << TColor::SavePrimitiveColor(GetFrameFillColor()) << ");\n";
2132 out << prefix << "SetFrameLineColor(" << TColor::SavePrimitiveColor(GetFrameLineColor()) << ");\n";
2133 out << prefix << "SetFrameFillStyle(" << GetFrameFillStyle() << ");\n";
2134 out << prefix << "SetFrameLineStyle(" << GetFrameLineStyle() << ");\n";
2135 out << prefix << "SetFrameLineWidth(" << GetFrameLineWidth() << ");\n";
2136 out << prefix << "SetFrameBorderSize(" << GetFrameBorderSize() << ");\n";
2137 out << prefix << "SetFrameBorderMode(" << GetFrameBorderMode() << ");\n";
2138 out << prefix << "SetHistFillColor(" << TColor::SavePrimitiveColor(GetHistFillColor()) << ");\n";
2139 out << prefix << "SetHistLineColor(" << TColor::SavePrimitiveColor(GetHistLineColor()) << ");\n";
2140 out << prefix << "SetHistFillStyle(" << GetHistFillStyle() << ");\n";
2141 out << prefix << "SetHistLineStyle(" << GetHistLineStyle() << ");\n";
2142 out << prefix << "SetHistLineWidth(" << GetHistLineWidth() << ");\n";
2143 out << prefix << "SetHistMinimumZero(" << asBool(GetHistMinimumZero()) << ");\n";
2144 out << prefix << "SetCanvasPreferGL(" << asBool(GetCanvasPreferGL()) << ");\n";
2145 out << prefix << "SetCanvasColor(" << TColor::SavePrimitiveColor(GetCanvasColor()) << ");\n";
2146 out << prefix << "SetCanvasBorderSize(" << GetCanvasBorderSize() << ");\n";
2147 out << prefix << "SetCanvasBorderMode(" << GetCanvasBorderMode() << ");\n";
2148 out << prefix << "SetCanvasDefH(" << GetCanvasDefH() << ");\n";
2149 out << prefix << "SetCanvasDefW(" << GetCanvasDefW() << ");\n";
2150 out << prefix << "SetCanvasDefX(" << GetCanvasDefX() << ");\n";
2151 out << prefix << "SetCanvasDefY(" << GetCanvasDefY() << ");\n";
2152 out << prefix << "SetPadColor(" << TColor::SavePrimitiveColor(GetPadColor()) << ");\n";
2153 out << prefix << "SetPadBorderSize(" << GetPadBorderSize() << ");\n";
2154 out << prefix << "SetPadBorderMode(" << GetPadBorderMode() << ");\n";
2155 out << prefix << "SetPadBottomMargin(" << GetPadBottomMargin() << ");\n";
2156 out << prefix << "SetPadTopMargin(" << GetPadTopMargin() << ");\n";
2157 out << prefix << "SetPadLeftMargin(" << GetPadLeftMargin() << ");\n";
2158 out << prefix << "SetPadRightMargin(" << GetPadRightMargin() << ");\n";
2159 out << prefix << "SetPadGridX(" << asBool(GetPadGridX()) << ");\n";
2160 out << prefix << "SetPadGridY(" << asBool(GetPadGridY()) << ");\n";
2161 out << prefix << "SetPadTickX(" << GetPadTickX() << ");\n";
2162 out << prefix << "SetPadTickY(" << GetPadTickY() << ");\n";
2163 out << prefix << "SetOrthoCamera(" << asBool(GetOrthoCamera()) << ");\n";
2164
2165 out << prefix << "SetCandleWhiskerRange(" << GetCandleWhiskerRange() << ");\n";
2166 out << prefix << "SetCandleBoxRange(" << GetCandleBoxRange() << ");\n";
2167 out << prefix << "SetCandleScaled(" << asBool(GetCandleScaled()) << ");\n";
2168 out << prefix << "SetViolinScaled(" << asBool(GetViolinScaled()) << ");\n";
2169
2170 // fPaperSizeX, fPaperSizeY
2171 out << prefix << "SetPaperSize(" << fPaperSizeX << ", " << fPaperSizeY << ");\n";
2172
2173 out << prefix << "SetScreenFactor(" << GetScreenFactor() << ");\n";
2174 out << prefix << "SetStatColor(" << TColor::SavePrimitiveColor(GetStatColor()) << ");\n";
2175 out << prefix << "SetStatTextColor(" << TColor::SavePrimitiveColor(GetStatTextColor()) << ");\n";
2176 out << prefix << "SetStatBorderSize(" << GetStatBorderSize() << ");\n";
2177 out << prefix << "SetStatFont(" << GetStatFont() << ");\n";
2178 out << prefix << "SetStatFontSize(" << GetStatFontSize() << ");\n";
2179 out << prefix << "SetStatStyle(" << GetStatStyle() << ");\n";
2180 out << prefix << "SetStatFormat(\"" << GetStatFormat() << "\");\n";
2181 out << prefix << "SetStatX(" << GetStatX() << ");\n";
2182 out << prefix << "SetStatY(" << GetStatY() << ");\n";
2183 out << prefix << "SetStatW(" << GetStatW() << ");\n";
2184 out << prefix << "SetStatH(" << GetStatH() << ");\n";
2185 out << prefix << "SetStripDecimals(" << asBool(GetStripDecimals()) << ");\n";
2186 out << prefix << "SetTitleAlign(" << GetTitleAlign() << ");\n";
2187 out << prefix << "SetTitleFillColor(" << TColor::SavePrimitiveColor(GetTitleFillColor()) << ");\n";
2188 out << prefix << "SetTitleTextColor(" << TColor::SavePrimitiveColor(GetTitleTextColor()) << ");\n";
2189 out << prefix << "SetTitleBorderSize(" << GetTitleBorderSize() << ");\n";
2190 out << prefix << "SetTitleFont(" << GetTitleFont() << ");\n";
2191 out << prefix << "SetTitleFontSize(" << GetTitleFontSize() << ");\n";
2192 out << prefix << "SetTitleStyle(" << GetTitleStyle() << ");\n";
2193 out << prefix << "SetTitleX(" << GetTitleX() << ");\n";
2194 out << prefix << "SetTitleY(" << GetTitleY() << ");\n";
2195 out << prefix << "SetTitleW(" << GetTitleW() << ");\n";
2196 out << prefix << "SetTitleH(" << GetTitleH() << ");\n";
2197 out << prefix << "SetLegoInnerR(" << GetLegoInnerR() << ");\n";
2198
2199 out << " \n";
2200
2201 // Palette belongs to the TColor
2202 // but stored together with any TStyle instance
2204
2205 out << " \n";
2206
2207 // fLineStyle
2208 for (Int_t li = 0; li < 30; ++li)
2209 out << prefix << "SetLineStyleString(" << li << ", \"" << TString(fLineStyle[li]).ReplaceSpecialCppChars()
2210 << "\");\n";
2211
2212 out << " \n";
2213
2214 out << prefix << "SetHeaderPS(\"" << TString(GetHeaderPS()).ReplaceSpecialCppChars() << "\");\n";
2215 out << prefix << "SetTitlePS(\"" << TString(GetTitlePS()).ReplaceSpecialCppChars() << "\");\n";
2216 out << prefix << "SetFitFormat(\"" << TString(GetFitFormat()).ReplaceSpecialCppChars() << "\");\n";
2217 out << prefix << "SetPaintTextFormat(\"" << TString(GetPaintTextFormat()).ReplaceSpecialCppChars() << "\");\n";
2218 out << prefix << "SetLineScalePS(" << GetLineScalePS() << ");\n";
2219 out << prefix << "SetJoinLinePS(" << GetJoinLinePS() << ");\n";
2220 out << prefix << "SetCapLinePS(" << GetCapLinePS() << ");\n";
2221 out << prefix << "SetColorModelPS(" << GetColorModelPS() << ");\n";
2222 out << prefix << "SetTimeOffset(" << TString::Format("%9.0f", GetTimeOffset()) << ");\n";
2223
2224 // Inheritance :
2225 SaveLineAttributes(out, "tmpStyle", -1, -1, -1);
2226 SaveFillAttributes(out, "tmpStyle", -1, -1);
2227 SaveMarkerAttributes(out, "tmpStyle", -1, -1, -1);
2228 SaveTextAttributes(out, "tmpStyle", 0, 0, 0, 0, 0);
2229}
#define c(i)
Definition RSha256.hxx:101
#define a(i)
Definition RSha256.hxx:99
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
short Style_t
Style number (short)
Definition RtypesCore.h:96
int Int_t
Signed integer 4 bytes (int)
Definition RtypesCore.h:59
short Color_t
Color number (short)
Definition RtypesCore.h:99
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
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
#define BIT(n)
Definition Rtypes.h:91
@ kRed
Definition Rtypes.h:67
@ kBlack
Definition Rtypes.h:66
@ kGreen
Definition Rtypes.h:67
@ kBlue
Definition Rtypes.h:67
R__EXTERN TApplication * gApplication
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
R__EXTERN TEnv * gEnv
Definition TEnv.h:170
Option_t Option_t option
Option_t Option_t SetLineWidth
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 SetTextSize
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 Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h offset
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 np
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 Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h length
Option_t Option_t SetTextFont
Option_t Option_t TPoint TPoint const char mode
Option_t Option_t SetFillColor
Option_t Option_t SetMarkerStyle
Option_t Option_t style
Option_t Option_t TPoint TPoint const char text
char name[80]
Definition TGX11.cxx:110
R__EXTERN TVirtualMutex * gROOTMutex
Definition TROOT.h:63
#define gROOT
Definition TROOT.h:411
TStyle * gStyle
Definition TStyle.cxx:30
const UInt_t kTakeStyle
Definition TStyle.cxx:31
R__EXTERN TStyle * gStyle
Definition TStyle.h:442
R__EXTERN TSystem * gSystem
Definition TSystem.h:572
#define R__LOCKGUARD(mutex)
#define gPad
Color * colors
Definition X3DBuffer.c:21
#define snprintf
Definition civetweb.c:1579
char ** Argv() const
virtual Color_t GetTitleColor() const
Definition TAttAxis.h:47
virtual Color_t GetLabelColor() const
Definition TAttAxis.h:39
virtual Int_t GetNdivisions() const
Definition TAttAxis.h:37
virtual Color_t GetAxisColor() const
Definition TAttAxis.h:38
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition TAttAxis.cxx:279
virtual Style_t GetTitleFont() const
Definition TAttAxis.h:48
virtual Float_t GetLabelOffset() const
Definition TAttAxis.h:41
virtual void SetAxisColor(Color_t color=1, Float_t alpha=1.)
Set color of the line axis and tick marks.
Definition TAttAxis.cxx:141
virtual void SetLabelSize(Float_t size=0.04)
Set size of axis labels.
Definition TAttAxis.cxx:184
virtual Style_t GetLabelFont() const
Definition TAttAxis.h:40
virtual void SetTitleFont(Style_t font=62)
Set the title font.
Definition TAttAxis.cxx:308
virtual void SetLabelOffset(Float_t offset=0.005)
Set distance between the axis and the labels.
Definition TAttAxis.cxx:172
virtual void SetLabelFont(Style_t font=62)
Set labels' font.
Definition TAttAxis.cxx:161
virtual void SetTitleSize(Float_t size=0.04)
Set size of axis title.
Definition TAttAxis.cxx:290
virtual void SetTitleColor(Color_t color=1)
Set color of axis title.
Definition TAttAxis.cxx:299
virtual Float_t GetTitleSize() const
Definition TAttAxis.h:45
virtual Float_t GetLabelSize() const
Definition TAttAxis.h:42
virtual Float_t GetTickLength() const
Definition TAttAxis.h:46
virtual void ResetAttAxis(Option_t *option="")
Reset axis attributes.
Definition TAttAxis.cxx:78
virtual Float_t GetTitleOffset() const
Definition TAttAxis.h:44
virtual void SetTickLength(Float_t length=0.03)
Set tick mark length.
Definition TAttAxis.cxx:265
virtual void SetNdivisions(Int_t n=510, Bool_t optim=kTRUE)
Set the number of divisions for this axis.
Definition TAttAxis.cxx:214
void Copy(TAttAxis &attaxis) const
Copy of the object.
Definition TAttAxis.cxx:60
virtual void SetLabelColor(Color_t color=1, Float_t alpha=1.)
Set color of labels.
Definition TAttAxis.cxx:151
Fill Area Attributes class.
Definition TAttFill.h:20
void Copy(TAttFill &attfill) const
Copy this fill attributes to a new TAttFill.
Definition TAttFill.cxx:206
virtual void ResetAttFill(Option_t *option="")
Reset this fill attributes to default values.
Definition TAttFill.cxx:229
virtual void SaveFillAttributes(std::ostream &out, const char *name, Int_t coldef=1, Int_t stydef=1001)
Save fill attributes as C++ statement(s) on output stream out.
Definition TAttFill.cxx:238
Line Attributes class.
Definition TAttLine.h:20
virtual void ResetAttLine(Option_t *option="")
Reset this line attributes to default values.
Definition TAttLine.cxx:264
void Copy(TAttLine &attline) const
Copy this line attributes to a new TAttLine.
Definition TAttLine.cxx:176
virtual void SaveLineAttributes(std::ostream &out, const char *name, Int_t coldef=1, Int_t stydef=1, Int_t widdef=1)
Save line attributes as C++ statement(s) on output stream out.
Definition TAttLine.cxx:274
Marker Attributes class.
Definition TAttMarker.h:20
virtual void SaveMarkerAttributes(std::ostream &out, const char *name, Int_t coldef=1, Int_t stydef=1, Int_t sizdef=1)
Save line attributes as C++ statement(s) on output stream out.
virtual void ResetAttMarker(Option_t *toption="")
Reset this marker attributes to the default values.
void Copy(TAttMarker &attmarker) const
Copy this marker attributes to a new TAttMarker.
virtual void SetMarkerSize(Size_t msize=1)
Set the marker size.
Definition TAttMarker.h:46
Text Attributes class.
Definition TAttText.h:20
virtual Float_t GetTextSize() const
Return the text size.
Definition TAttText.h:38
virtual void SetTextAlign(Short_t align=11)
Set the text alignment.
Definition TAttText.h:44
virtual Short_t GetTextAlign() const
Return the text alignment.
Definition TAttText.h:34
virtual Font_t GetTextFont() const
Return the text font.
Definition TAttText.h:37
virtual Color_t GetTextColor() const
Return the text color.
Definition TAttText.h:36
virtual Float_t GetTextAngle() const
Return the text angle.
Definition TAttText.h:35
virtual void SaveTextAttributes(std::ostream &out, const char *name, Int_t alidef=12, Float_t angdef=0, Int_t coldef=1, Int_t fondef=61, Float_t sizdef=1)
Save text attributes as C++ statement(s) on output stream out.
Definition TAttText.cxx:372
virtual void SetTextSize(Float_t tsize=1)
Set the text size.
Definition TAttText.h:49
void Copy(TAttText &atttext) const
Copy this text attributes to a new TAttText.
Definition TAttText.cxx:293
virtual void ResetAttText(Option_t *toption="")
Reset this text attributes to default values.
Definition TAttText.cxx:360
Using a TBrowser one can browse all ROOT objects.
Definition TBrowser.h:37
The color creation and management class.
Definition TColor.h:22
static void SetPalette(Int_t ncolors, Int_t *colors, Float_t alpha=1.)
Static function.
Definition TColor.cxx:2933
static Int_t GetColorPalette(Int_t i)
Static function returning the color number i in current palette.
Definition TColor.cxx:1509
static void SaveColorsPalette(std::ostream &out)
Store current palette in the output macro.
Definition TColor.cxx:3648
static Int_t CreateColorTableFromFile(TString fileName, Float_t alpha=1.)
Static function creating a color palette based on an input text file.
Definition TColor.cxx:2636
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 Int_t GetNumberOfColors()
Static function returning number of colors in the color palette.
Definition TColor.cxx:1529
This class stores the date and time with a precision of one second in an unsigned 32 bit word (950130...
Definition TDatime.h:37
const char * AsSQLString() const
Return the date & time in SQL compatible string format, like: 1997-01-15 20:16:28.
Definition TDatime.cxx:151
virtual Int_t GetValue(const char *name, Int_t dflt) const
Returns the integer value for a resource.
Definition TEnv.cxx:490
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
const char * GetTitle() const override
Returns title of object.
Definition TNamed.h:50
virtual void SetNameTitle(const char *name, const char *title)
Set all the TNamed parameters (name and title).
Definition TNamed.cxx:163
Mother of all ROOT objects.
Definition TObject.h:41
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:1071
static void SavePrimitiveConstructor(std::ostream &out, TClass *cl, const char *variable_name, const char *constructor_agrs="", Bool_t empty_line=kTRUE)
Save object constructor in the output stream "out".
Definition TObject.cxx:771
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
const char * Data() const
Definition TString.h:384
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
TStyle objects may be created to define special styles.
Definition TStyle.h:29
Bool_t GetViolinScaled() const
Definition TStyle.h:295
Double_t GetTimeOffset() const
Definition TStyle.h:271
Int_t GetAxisMaxDigits() const
Method returns maximum number of digits permitted for the axis labels above which the notation with 1...
Definition TStyle.cxx:1158
Int_t GetOptLogy() const
Definition TStyle.h:250
Color_t GetGridColor() const
Definition TStyle.h:224
Int_t fOptFile
True if option File is selected.
Definition TStyle.h:45
Int_t GetOptStat() const
Definition TStyle.h:247
Color_t GetLabelColor(Option_t *axis="X") const
Return the label color number in the axis.
Definition TStyle.cxx:1110
void SetAxisColor(Color_t color=1, Option_t *axis="X")
Set color to draw the axis line and tick marks.
Definition TStyle.cxx:1324
void SetPadBorderMode(Int_t mode=1)
Definition TStyle.h:361
void SetOptTitle(Int_t tit=1)
Definition TStyle.h:338
Color_t fGridColor
Grid line color (if 0 use axis line color)
Definition TStyle.h:60
void SetPadTopMargin(Float_t margin=0.1)
Definition TStyle.h:363
Color_t GetStatTextColor() const
Definition TStyle.h:260
void SetLegendFont(Style_t font=62)
Definition TStyle.h:357
void SetTitleX(Float_t x=0)
Definition TStyle.h:417
void SetOptStat(Int_t stat=1)
The type of information printed in the histogram statistics box can be selected via the parameter mod...
Definition TStyle.cxx:1641
Int_t fCanvasDefW
Default canvas width.
Definition TStyle.h:89
Int_t fHatchesLineWidth
Hatches line width for hatch styles > 3100.
Definition TStyle.h:68
Float_t fScreenFactor
Multiplication factor for canvas size and position.
Definition TStyle.h:105
Float_t fYAxisExpXOffset
Y axis exponent label X offset.
Definition TStyle.h:149
Float_t fTitleFontSize
Font size in pixels for fonts with precision type 3.
Definition TStyle.h:123
void SetDateX(Float_t x=0.01)
Definition TStyle.h:341
Int_t fCanvasDefX
Default canvas top X position.
Definition TStyle.h:90
void SetStatFormat(const char *format="6.4g")
Definition TStyle.h:400
void SetPadBottomMargin(Float_t margin=0.1)
Definition TStyle.h:362
Float_t GetTitleX() const
Definition TStyle.h:282
Bool_t fViolinScaled
Violin plot, shall the violin or histos be scaled to each other by the maximum height?
Definition TStyle.h:144
void SetLegendFillColor(Color_t color=0)
Definition TStyle.h:355
Float_t fPadRightMargin
Pad right margin.
Definition TStyle.h:98
Float_t fPaperSizeX
PostScript paper size along X.
Definition TStyle.h:103
void SetCandleCrossLineWidth(Int_t CrossLineWidth=1)
Set the line width of the cross marker of a candle plot ([1,5]).
Definition TStyle.cxx:1983
void SetPaintTextFormat(const char *format="g")
Definition TStyle.h:390
Int_t GetOptTitle() const
Definition TStyle.h:248
Float_t GetScreenFactor() const
Definition TStyle.h:258
TString fPaintTextFormat
Printing format for TH2::PaintText.
Definition TStyle.h:134
Color_t GetHistLineColor() const
Definition TStyle.h:235
Color_t fStatColor
Stat fill area color.
Definition TStyle.h:106
Int_t GetNdivisions(Option_t *axis="X") const
Return number of divisions.
Definition TStyle.cxx:1078
void SetFrameFillColor(Color_t color=1)
Definition TStyle.h:376
Int_t fFrameBorderMode
Pad frame border mode.
Definition TStyle.h:76
Int_t GetPadTickX() const
Definition TStyle.h:219
Color_t GetTitleColor(Option_t *axis="X") const
Return title color.
Definition TStyle.cxx:1205
Style_t fHistLineStyle
Histogram line style.
Definition TStyle.h:80
Color_t GetFrameLineColor() const
Definition TStyle.h:228
Style_t GetGridStyle() const
Definition TStyle.h:225
Int_t GetJoinLinePS() const
Returns the line join method used for PostScript, PDF and SVG output. See TPostScript::SetLineJoin fo...
Definition TStyle.h:289
Style_t fLegendFillStyle
Legend fill style.
Definition TStyle.h:65
void SetStatFont(Style_t font=62)
Definition TStyle.h:398
void SetEndErrorSize(Float_t np=2)
Set the size (in pixels) of the small lines drawn at the end of the error bars (TH1 or TGraphErrors).
Definition TStyle.cxx:1341
Float_t GetStatFontSize() const
Definition TStyle.h:263
Float_t GetBarOffset() const
Definition TStyle.h:184
Float_t GetStatX() const
Definition TStyle.h:266
Float_t GetLabelSize(Option_t *axis="X") const
Return label size.
Definition TStyle.cxx:1146
void SetCapLinePS(Int_t capline=0)
Set the line cap method used for PostScript, PDF and SVG output. See TPostScript::SetLineCap for deta...
Definition TStyle.h:311
Color_t GetPadColor() const
Definition TStyle.h:210
Float_t fBarOffset
Offset of bar for graphs.
Definition TStyle.h:36
Bool_t GetCandleScaled() const
Definition TStyle.h:294
void SetPadRightMargin(Float_t margin=0.1)
Definition TStyle.h:365
Style_t GetHistFillStyle() const
Definition TStyle.h:236
void SetCanvasColor(Color_t color=19)
Definition TStyle.h:347
TAttAxis fZaxis
Z axis attributes.
Definition TStyle.h:34
void SetTitleFont(Style_t font=62, Option_t *axis="X")
Definition TStyle.cxx:1776
Float_t GetPadRightMargin() const
Definition TStyle.h:216
Float_t GetTickLength(Option_t *axis="X") const
Return tick length.
Definition TStyle.cxx:1193
void SetExponentOffset(Float_t xoff=0., Float_t yoff=0., Option_t *axis="XY")
Method set X and Y offset of the axis 10^n notation.
Definition TStyle.cxx:1837
Int_t fDrawBorder
Flag to draw border(=1) or not (0)
Definition TStyle.h:38
void SetTitleBorderSize(Width_t size=2)
Definition TStyle.h:412
void SaveSource(const char *filename, Option_t *option=nullptr)
Save the current style in a C++ macro file.
Definition TStyle.cxx:1997
Double_t GetCandleWhiskerRange() const
Definition TStyle.h:292
Style_t GetFrameFillStyle() const
Definition TStyle.h:229
Float_t GetTitleSize(Option_t *axis="X") const
Return title size.
Definition TStyle.cxx:1241
Float_t GetLegoInnerR() const
Definition TStyle.h:242
Int_t AxisChoice(Option_t *axis) const
Definition TStyle.cxx:506
TString fTitlePS
User defined Postscript file title.
Definition TStyle.h:132
Style_t GetLabelFont(Option_t *axis="X") const
Return label font.
Definition TStyle.cxx:1122
virtual ~TStyle()
Destructor.
Definition TStyle.cxx:478
void SetCanvasBorderMode(Int_t mode=1)
Definition TStyle.h:349
Float_t GetTitleY() const
Definition TStyle.h:283
Style_t fFuncStyle
Function style.
Definition TStyle.h:58
void SetOptDate(Int_t datefl=1)
If optdate is non null, the current date/time will be printed in the canvas.
Definition TStyle.cxx:1547
Int_t fPadTickX
True to set special pad ticks along X.
Definition TStyle.h:101
Float_t GetDateX() const
Definition TStyle.h:199
Style_t fStatFont
Font style of Stats PaveLabel.
Definition TStyle.h:109
TString fLineStyle[30]
String describing line style i (for postScript)
Definition TStyle.h:130
void SetDateY(Float_t y=0.01)
Definition TStyle.h:342
Width_t fGridWidth
Grid line width.
Definition TStyle.h:62
Double_t fCandleBoxRange
Candle plot, The fraction which is covered by the box (0 < x < 1), default 0.5.
Definition TStyle.h:142
Float_t GetTitleOffset(Option_t *axis="X") const
Return title offset.
Definition TStyle.cxx:1229
Color_t GetHistFillColor() const
Definition TStyle.h:234
void SetFrameBorderMode(Int_t mode=1)
Definition TStyle.h:382
Style_t GetTitleFont(Option_t *axis="X") const
Return title font.
Definition TStyle.cxx:1217
Int_t fOptLogx
True if log scale in X.
Definition TStyle.h:39
Bool_t fStripDecimals
Strip decimals in axis labels.
Definition TStyle.h:117
Bool_t fHistMinimumZero
True if default minimum is 0, false if minimum is automatic.
Definition TStyle.h:82
Style_t fLegendFont
Legend font style.
Definition TStyle.h:66
void SetFuncColor(Color_t color=1)
Definition TStyle.h:371
Bool_t GetHistMinimumZero() const
Definition TStyle.h:239
TAttAxis fXaxis
X axis attributes.
Definition TStyle.h:32
void SetHeaderPS(const char *header)
Define a string to be inserted in the Postscript header.
Definition TStyle.cxx:1363
Float_t GetStatY() const
Definition TStyle.h:267
Color_t fTitleColor
Title fill area color.
Definition TStyle.h:119
Int_t fJoinLinePS
Determines the appearance of joining lines on PostScript, PDF and SVG.
Definition TStyle.h:136
void SetPadTickY(Int_t ticky)
Definition TStyle.h:369
void SetTitleOffset(Float_t offset=1, Option_t *axis="X")
Specify a parameter offset to control the distance between the axis and the axis title.
Definition TStyle.cxx:1799
Width_t fFrameLineWidth
Pad frame line width.
Definition TStyle.h:74
Color_t GetTitleFillColor() const
Definition TStyle.h:273
Int_t GetCanvasDefH() const
Definition TStyle.h:193
void SetColorModelPS(Int_t c=0)
Define the color model used by TPostScript and TPDF (RGB or CMYK).
Definition TStyle.cxx:1276
Style_t GetTitleStyle() const
Definition TStyle.h:275
void SetPadTickX(Int_t tickx)
Definition TStyle.h:368
Float_t GetLabelOffset(Option_t *axis="X") const
Return label offset.
Definition TStyle.cxx:1134
Int_t GetCanvasDefX() const
Definition TStyle.h:195
Int_t fOptTitle
True if option Title is selected.
Definition TStyle.h:44
Style_t fTitleFont
Font style of Title PaveLabel.
Definition TStyle.h:122
Int_t GetOptDate() const
Definition TStyle.h:244
Bool_t GetPadGridY() const
Definition TStyle.h:218
Int_t fOptDate
True if date option is selected.
Definition TStyle.h:42
Float_t fPadTopMargin
Pad top margin.
Definition TStyle.h:96
Color_t GetStatColor() const
Definition TStyle.h:259
Width_t fCanvasBorderSize
Canvas border size.
Definition TStyle.h:86
void SetTitleTextColor(Color_t color=1)
Definition TStyle.h:409
void SetAxisMaxDigits(Int_t maxd=5)
Method set maximum number of digits permitted for the axis labels above which the notation with 10^N ...
Definition TStyle.cxx:1881
Color_t fCanvasColor
Canvas color.
Definition TStyle.h:85
Float_t GetPadLeftMargin() const
Definition TStyle.h:215
Int_t fCandleCrossLineWidth
Line width of the cross marker of a candle plot ([1,5]).
Definition TStyle.h:146
Double_t GetHatchesSpacing() const
Definition TStyle.h:203
void SetPalette(Int_t ncolors=kBird, Int_t *colors=nullptr, Float_t alpha=1.)
See TColor::SetPalette.
Definition TStyle.cxx:1889
Float_t fStatY
Y position of top right corner of stat box.
Definition TStyle.h:114
void Copy(TObject &style) const override
Copy this style.
Definition TStyle.cxx:550
Width_t GetLegendBorderSize() const
Definition TStyle.h:204
Style_t fFrameLineStyle
Pad frame line style.
Definition TStyle.h:73
void SetStatBorderSize(Width_t size=2)
Definition TStyle.h:397
Float_t GetBarWidth() const
Definition TStyle.h:185
Int_t fCanvasBorderMode
Canvas border mode.
Definition TStyle.h:87
Bool_t GetCanvasPreferGL() const
Definition TStyle.h:189
Int_t GetColorModelPS() const
Definition TStyle.h:198
void SetErrorX(Float_t errorx=0.5)
Definition TStyle.h:344
Int_t GetCanvasDefY() const
Definition TStyle.h:196
void SetTitleColor(Color_t color=1, Option_t *axis="X")
Definition TStyle.cxx:1755
void SetNumberContours(Int_t number=20)
Set the default number of contour levels when drawing 2-d plots.
Definition TStyle.cxx:1501
void Paint(Option_t *option="") override
Copy this style to gStyle.
Definition TStyle.cxx:1253
Float_t fTitleX
X position of top left corner of title box.
Definition TStyle.h:125
Int_t fShowToolBar
Show toolbar.
Definition TStyle.h:49
void SetLabelFont(Style_t font=62, Option_t *axis="X")
Set font number used to draw axis labels.
Definition TStyle.cxx:1414
Int_t fCapLinePS
Determines the appearance of line caps on PostScript, PDF and SVG.
Definition TStyle.h:137
Width_t fFuncWidth
Function line width.
Definition TStyle.h:59
Int_t fPadTickY
True to set special pad ticks along Y.
Definition TStyle.h:102
TAttAxis fYaxis
Y axis attributes.
Definition TStyle.h:33
Float_t fImageScaling
Image scaling to produce high definition bitmap images.
Definition TStyle.h:140
Double_t fCandleWhiskerRange
Candle plot, the fraction which is covered by the whiskers (0 < x < 1), default 1.
Definition TStyle.h:141
Width_t GetFrameBorderSize() const
Definition TStyle.h:232
Width_t fHistLineWidth
Histogram line width.
Definition TStyle.h:81
void SetTimeOffset(Double_t toffset)
Change the time offset for time plotting.
Definition TStyle.cxx:1914
Bool_t fPadGridY
True to get the grid along Y.
Definition TStyle.h:100
void SetTitlePS(const char *pstitle)
Define a string to be used in the %Title of the Postscript files.
Definition TStyle.cxx:1384
Style_t fHistFillStyle
Histogram fill style.
Definition TStyle.h:79
void SetHistMinimumZero(Bool_t zero=kTRUE)
If the argument zero=kTRUE the minimum value for the Y axis of 1-d histograms is set to 0.
Definition TStyle.cxx:1289
Bool_t GetPadGridX() const
Definition TStyle.h:217
Float_t GetStatH() const
Definition TStyle.h:269
void SetPadLeftMargin(Float_t margin=0.1)
Definition TStyle.h:364
Float_t fEndErrorSize
Size of lines at the end of error bars.
Definition TStyle.h:55
void SetJoinLinePS(Int_t joinline=0)
Set the line join method used for PostScript, PDF and SVG output. See TPostScript::SetLineJoin for de...
Definition TStyle.h:310
static TClass * Class()
Bool_t fIsReading
! Set to FALSE when userclass::UseCurrentStyle is called by the style manager
Definition TStyle.h:139
Width_t GetGridWidth() const
Definition TStyle.h:226
Width_t fFrameBorderSize
Pad frame border size.
Definition TStyle.h:75
Color_t GetFuncColor() const
Definition TStyle.h:221
void SetTitleXOffset(Float_t offset=1)
Definition TStyle.h:413
void SetLegendBorderSize(Width_t size=4)
Definition TStyle.h:354
TAttText * GetAttDate()
Definition TStyle.h:170
Color_t fTitleTextColor
Title text color.
Definition TStyle.h:120
Int_t GetPadTickY() const
Definition TStyle.h:220
Width_t GetPadBorderSize() const
Definition TStyle.h:211
void SetStripDecimals(Bool_t strip=kTRUE)
Set option to strip decimals when drawing axis labels.
Definition TStyle.cxx:1927
Width_t GetTitleBorderSize() const
Definition TStyle.h:277
void SetHistLineColor(Color_t color=1)
Definition TStyle.h:384
TString fHeaderPS
User defined additional Postscript header.
Definition TStyle.h:131
Int_t GetColorPalette(Int_t i) const
Return color number i in current palette.
Definition TStyle.cxx:1102
virtual void cd()
Change current style.
Definition TStyle.cxx:542
const char * GetLineStyleString(Int_t i=1) const
Return line style string (used by PostScript).
Definition TStyle.cxx:1167
void SetLabelOffset(Float_t offset=0.005, Option_t *axis="X")
Set offset between axis and axis labels.
Definition TStyle.cxx:1430
Color_t fPadColor
Pad color.
Definition TStyle.h:92
void SetFitFormat(const char *format="5.4g")
Definition TStyle.h:305
Float_t GetErrorX() const
Definition TStyle.h:188
Float_t fDateX
X position of the date in the canvas (in NDC)
Definition TStyle.h:53
Int_t fOptLogz
True if log scale in z.
Definition TStyle.h:41
void SetTitleSize(Float_t size=0.02, Option_t *axis="X")
Definition TStyle.cxx:1818
void SetTitleFillColor(Color_t color=1)
Definition TStyle.h:408
TString fFitFormat
Printing format for fit parameters.
Definition TStyle.h:133
Int_t fPadBorderMode
Pad border mode.
Definition TStyle.h:94
Double_t GetLegendTextSize() const
Definition TStyle.h:208
Int_t fNumberContours
Default number of contours for 2-d plots.
Definition TStyle.h:51
void SetLineStyleString(Int_t i, const char *text)
Set line style string using the PostScript convention.
Definition TStyle.cxx:1484
Double_t GetCandleBoxRange() const
Definition TStyle.h:293
Float_t fPadLeftMargin
Pad left margin.
Definition TStyle.h:97
Float_t fTitleY
Y position of top left corner of title box.
Definition TStyle.h:126
Double_t fTimeOffset
Time offset to the beginning of an axis.
Definition TStyle.h:138
Color_t GetCanvasColor() const
Definition TStyle.h:190
Color_t fFrameLineColor
Pad frame line color.
Definition TStyle.h:71
Double_t fLegendTextSize
Legend text size. If 0 the size is computed automatically.
Definition TStyle.h:67
Int_t fShowEditor
Show pad editor.
Definition TStyle.h:48
Color_t fLegendFillColor
Legend fill color.
Definition TStyle.h:64
void SetTitleAlign(Int_t a=13)
Definition TStyle.h:407
Float_t fLineScalePS
Line scale factor when drawing lines on Postscript.
Definition TStyle.h:135
void SetPaperSize(EPaperSize size)
Set paper size for PostScript output.
Definition TStyle.cxx:1705
Float_t GetEndErrorSize() const
Definition TStyle.h:187
Float_t GetPadBottomMargin() const
Definition TStyle.h:213
void SetFrameLineWidth(Width_t width=1)
Definition TStyle.h:380
void SetTickLength(Float_t length=0.03, Option_t *axis="X")
Set the tick marks length for an axis.
Definition TStyle.cxx:1736
Double_t fHistTopMargin
Margin between histogram's top and pad's top.
Definition TStyle.h:83
void SetNdivisions(Int_t n=510, Option_t *axis="X")
Set the number of divisions to draw an axis.
Definition TStyle.cxx:1310
Int_t fOptStat
True if option Stat is selected.
Definition TStyle.h:43
Double_t fHatchesSpacing
Hatches spacing for hatch styles > 3100.
Definition TStyle.h:69
Color_t fStatTextColor
Stat text color.
Definition TStyle.h:107
Width_t GetFuncWidth() const
Definition TStyle.h:223
TAttText fAttDate
Canvas date attribute.
Definition TStyle.h:52
Width_t fStatBorderSize
Border size of Stats PaveLabel.
Definition TStyle.h:108
Bool_t GetOrthoCamera() const
Definition TStyle.h:298
void Browse(TBrowser *b) override
Browse the style object.
Definition TStyle.cxx:516
Float_t fTitleW
Width of title box.
Definition TStyle.h:127
Width_t fPadBorderSize
Pad border size.
Definition TStyle.h:93
void SetFuncWidth(Width_t width=4)
Definition TStyle.h:372
void SetLegendFillStyle(Style_t style=1001)
Definition TStyle.h:356
Int_t GetDrawBorder() const
Definition TStyle.h:186
const char * GetTitlePS() const
Definition TStyle.h:287
Color_t fHistLineColor
Histogram line color.
Definition TStyle.h:78
Int_t fOptFit
True if option Fit is selected.
Definition TStyle.h:46
Int_t GetCanvasDefW() const
Definition TStyle.h:194
Int_t GetCapLinePS() const
Returns the line cap method used for PostScript, PDF and SVG output. See TPostScript::SetLineCap for ...
Definition TStyle.h:290
Int_t DistancetoPrimitive(Int_t px, Int_t py) override
Function used by the TStyle manager when drawing a canvas showing the current style.
Definition TStyle.cxx:686
Width_t GetStatBorderSize() const
Definition TStyle.h:261
void GetPaperSize(Float_t &xsize, Float_t &ysize) const
Set paper size for PostScript output.
Definition TStyle.cxx:1184
Int_t GetStripDecimals() const
Definition TStyle.h:270
Float_t fXAxisExpYOffset
X axis exponent label Y offset.
Definition TStyle.h:148
Style_t GetLegendFillStyle() const
Definition TStyle.h:206
Color_t fFrameFillColor
Pad frame fill color.
Definition TStyle.h:70
Style_t GetHistLineStyle() const
Definition TStyle.h:237
TStyle()
Default constructor.
Definition TStyle.cxx:144
Float_t fTitleH
Height of title box.
Definition TStyle.h:128
void SetTitleStyle(Style_t style=1001)
Definition TStyle.h:410
void SetStatColor(Color_t color=19)
Definition TStyle.h:394
Style_t fFrameFillStyle
Pad frame fill style.
Definition TStyle.h:72
void SetPadColor(Color_t color=19)
Definition TStyle.h:359
virtual void Reset(Option_t *option="")
Reset.
Definition TStyle.cxx:695
Int_t fShowEventStatus
Show event status panel.
Definition TStyle.h:47
Color_t GetTitleTextColor() const
Definition TStyle.h:274
Bool_t fPadGridX
True to get the grid along X.
Definition TStyle.h:99
void SetStatY(Float_t y=0)
Definition TStyle.h:402
void SetCandleWhiskerRange(Double_t wRange=1.0)
By setting whisker-range for candle plot, one can force the whiskers to cover the fraction of the dis...
Definition TStyle.cxx:1940
Style_t fTitleStyle
Fill area style of title PaveLabel.
Definition TStyle.h:124
Style_t GetLegendFont() const
Definition TStyle.h:207
Int_t GetOptLogx() const
Definition TStyle.h:249
Int_t fCandleCircleLineWidth
Line width of the circle marker of a candle plot ([1,5]).
Definition TStyle.h:145
void SetLegendTextSize(Double_t size=0.)
Definition TStyle.h:358
TStyle & operator=(const TStyle &style)
Assignment operator.
Definition TStyle.cxx:496
void SavePrimitive(std::ostream &out, Option_t *="") override
Save primitive as a C++ statement(s) on output stream out.
Definition TStyle.cxx:2042
Bool_t fCandleScaled
Candle plot, shall the box-width be scaled to each other by the integral of a box?
Definition TStyle.h:143
Float_t GetTitleH() const
Definition TStyle.h:285
Style_t GetStatStyle() const
Definition TStyle.h:264
Float_t fStatFontSize
Font size in pixels for fonts with precision type 3.
Definition TStyle.h:110
Bool_t fCanvasPreferGL
If true, rendering in canvas is with GL.
Definition TStyle.h:84
Float_t fLegoInnerR
Inner radius for cylindrical legos.
Definition TStyle.h:129
Width_t GetHistLineWidth() const
Definition TStyle.h:238
EPaperSize
Definition TStyle.h:155
@ kA4
Definition TStyle.h:155
@ kUSLetter
Definition TStyle.h:155
Int_t fOptLogy
True if log scale in y.
Definition TStyle.h:40
Int_t fAxisMaxDigits
Number of digits above which the 10^N notation is used for axis.
Definition TStyle.h:151
Style_t GetFrameLineStyle() const
Definition TStyle.h:230
void SetIsReading(Bool_t reading=kTRUE)
Sets the fIsReading member to reading (default=kTRUE).
Definition TStyle.cxx:1375
Float_t GetStatW() const
Definition TStyle.h:268
Style_t fGridStyle
Grid line style.
Definition TStyle.h:61
Color_t fHistFillColor
Histogram fill color.
Definition TStyle.h:77
Float_t GetDateY() const
Definition TStyle.h:200
const char * GetFitFormat() const
Definition TStyle.h:201
Int_t fCanvasDefH
Default canvas height.
Definition TStyle.h:88
Int_t GetCanvasBorderMode() const
Definition TStyle.h:192
Int_t GetPadBorderMode() const
Definition TStyle.h:212
Float_t fXAxisExpXOffset
X axis exponent label X offset.
Definition TStyle.h:147
Color_t fFuncColor
Function color.
Definition TStyle.h:57
const char * GetHeaderPS() const
Definition TStyle.h:286
void SetTitleYOffset(Float_t offset=1)
Definition TStyle.h:415
const char * GetStatFormat() const
Definition TStyle.h:265
Width_t GetCanvasBorderSize() const
Definition TStyle.h:191
Int_t GetNumberOfColors() const
Return number of colors in the color palette.
Definition TStyle.cxx:1176
Float_t fStatH
Height of stat box.
Definition TStyle.h:116
Int_t GetOptFit() const
Definition TStyle.h:246
Float_t fStatX
X position of top right corner of stat box.
Definition TStyle.h:113
Int_t fCanvasDefY
Default canvas top Y position.
Definition TStyle.h:91
Int_t GetNumberContours() const
Definition TStyle.h:243
void SetHistLineWidth(Width_t width=1)
Definition TStyle.h:387
Bool_t fOrthoCamera
Use orthographic camera with web display.
Definition TStyle.h:152
Float_t fPadBottomMargin
Pad bottom margin.
Definition TStyle.h:95
const char * GetPaintTextFormat() const
Definition TStyle.h:252
Float_t GetLineScalePS() const
Definition TStyle.h:291
Float_t fErrorX
Per cent of bin width for errors along X.
Definition TStyle.h:56
void SetLabelColor(Color_t color=1, Option_t *axis="X")
Set axis labels color.
Definition TStyle.cxx:1394
void SetCandleBoxRange(Double_t bRange=0.5)
By setting box-range for candle plot, one can force the box of the candle-chart to cover that given f...
Definition TStyle.cxx:1957
Style_t GetStatFont() const
Definition TStyle.h:262
void SetLabelSize(Float_t size=0.04, Option_t *axis="X")
Set size of axis labels.
Definition TStyle.cxx:1445
void SetOptFit(Int_t fit=1)
The type of information about fit parameters printed in the histogram statistics box can be selected ...
Definition TStyle.cxx:1594
void SetCandleCircleLineWidth(Int_t CircleLineWidth=1)
Set the line width of the circle marker of a candle plot ([1,5]).
Definition TStyle.cxx:1970
Width_t fLegendBorderSize
Legend box border size.
Definition TStyle.h:63
Float_t fStatW
Width of stat box.
Definition TStyle.h:115
Float_t fYAxisExpYOffset
Y axis exponent label Y offset.
Definition TStyle.h:150
Float_t fBarWidth
Width of bar for graphs.
Definition TStyle.h:35
Width_t fTitleBorderSize
Border size of Title PavelLabel.
Definition TStyle.h:121
Int_t GetOptLogz() const
Definition TStyle.h:251
Float_t fPaperSizeY
PostScript paper size along Y.
Definition TStyle.h:104
static void BuildStyles()
Create some standard styles.
Definition TStyle.cxx:524
Int_t fTitleAlign
Title box alignment.
Definition TStyle.h:118
Style_t GetFuncStyle() const
Definition TStyle.h:222
Color_t GetLegendFillColor() const
Definition TStyle.h:205
Float_t GetTitleFontSize() const
Definition TStyle.h:276
Int_t GetHatchesLineWidth() const
Definition TStyle.h:202
Float_t fDateY
Y position of the date in the canvas (in NDC)
Definition TStyle.h:54
void GetExponentOffset(Float_t &xoff, Float_t &yoff, Option_t *axis="X") const
Method returns X and Y offset of the axis 10^n notation.
Definition TStyle.cxx:1856
Int_t GetTitleAlign() const
Definition TStyle.h:272
Color_t GetAxisColor(Option_t *axis="X") const
Return the axis color number in the axis.
Definition TStyle.cxx:1090
Int_t GetFrameBorderMode() const
Definition TStyle.h:233
TString fStatFormat
Printing format for stats.
Definition TStyle.h:112
Float_t GetPadTopMargin() const
Definition TStyle.h:214
Int_t fColorModelPS
PostScript color model: 0 = RGB, 1 = CMYK.
Definition TStyle.h:37
void SetLineScalePS(Float_t scale=3)
Definition TStyle.h:312
Width_t GetFrameLineWidth() const
Definition TStyle.h:231
Color_t GetFrameFillColor() const
Definition TStyle.h:227
Style_t fStatStyle
Fill area style of Stats PaveLabel.
Definition TStyle.h:111
Float_t GetTitleW() const
Definition TStyle.h:284
virtual const char * BaseName(const char *pathname)
Base name of a file name. Base name of /user/root is root.
Definition TSystem.cxx:944
const Int_t n
Definition legend1.C:16
TLine l
Definition textangle.C:4