Logo ROOT  
Reference Guide

The histogram painter class.

Implements all histograms' drawing's options.

  • Introduction
  • Histograms' plotting options
    • Options supported for 1D and 2D histograms
    • Options supported for 1D histograms
    • Options supported for 2D histograms
    • Options supported for 3D histograms
    • Options supported for histograms' stacks (THStack)
  • Setting the Style
  • Setting line, fill, marker, and text attributes
  • Setting Tick marks on the histogram axis
  • Giving titles to the X, Y and Z axis
  • The option "SAME"
    • Limitations
  • Colors automatically picked in palette
  • Superimposing two histograms with different scales in the same pad
  • Statistics Display
  • Fit Statistics
  • The error bars options
  • The bar chart option
  • The "BAR" and "HBAR" options
  • The SCATter plot option (default for 2D histograms)
  • The ARRow option
  • The BOX option
  • The COLor option
  • The CANDLE and VIOLIN options
    • The CANDLE option
    • The VIOLIN option
  • The TEXT and TEXTnn Option
  • The CONTour options
    • The LIST option
    • The AITOFF, MERCATOR, SINUSOIDAL and PARABOLIC options
  • The LEGO options
  • The "SURFace" options
  • Cylindrical, Polar, Spherical and PseudoRapidity/Phi options
  • Base line for bar-charts and lego plots
  • TH2Poly Drawing
  • The SPEC option
  • Option "Z" : Adding the color palette on the right side of the pad
  • Setting the color palette
  • Drawing a sub-range of a 2-D histogram; the [cutg] option
  • Drawing options for 3D histograms
  • Drawing option for histograms' stacks
  • Drawing of 3D implicit functions
  • Associated functions drawing
  • Drawing using OpenGL
    • General information: plot types and supported options
    • TH3 as color boxes
    • TH3 as boxes (spheres)
    • TH3 as iso-surface(s)
    • TF3 (implicit function)
    • Parametric surfaces
    • Interaction with the plots
    • Selectable parts
    • Rotation and zooming
    • Panning
    • Box cut
    • Plot specific interactions (dynamic slicing etc.)
    • Surface with option "GLSURF"
    • TF3
    • Box
    • Iso
    • Parametric plot
  • Highlight mode for histogram
    • Highlight mode and user function

Introduction

Histograms are drawn via the THistPainter class. Each histogram has a pointer to its own painter (to be usable in a multithreaded program). When the canvas has to be redrawn, the Paint function of each objects in the pad is called. In case of histograms, TH1::Paint invokes directly THistPainter::Paint.

To draw a histogram h it is enough to do:

h->Draw();

h can be of any kind: 1D, 2D or 3D. To choose how the histogram will be drawn, the Draw() method can be invoked with an option. For instance to draw a 2D histogram as a lego plot it is enough to do:

h->Draw("lego");

THistPainter offers many options to paint 1D, 2D and 3D histograms.

When the Draw() method of a histogram is called for the first time (TH1::Draw), it creates a THistPainter object and saves a pointer to this "painter" as a data member of the histogram. The THistPainter class specializes in the drawing of histograms. It is separated from the histogram so that one can have histograms without the graphics overhead, for example in a batch program. Each histogram having its own painter (rather than a central singleton painter painting all histograms), allows two histograms to be drawn in two threads without overwriting the painter's values.

When a displayed histogram is filled again, there is no need to call the Draw() method again; the image will be refreshed the next time the pad will be updated.

A pad is updated after one of these three actions:

  1. a carriage control on the ROOT command line,
  2. a click inside the pad,
  3. a call to TPad::Update.

By default a call to TH1::Draw() clears the pad of all objects before drawing the new image of the histogram. One can use the SAME option to leave the previous display intact and superimpose the new histogram. The same histogram can be drawn with different graphics options in different pads.

When a displayed histogram is deleted, its image is automatically removed from the pad.

To create a copy of the histogram when drawing it, one can use TH1::DrawClone(). This will clone the histogram and allow to change and delete the original one without affecting the clone.

Histograms' plotting options

Most options can be concatenated with or without spaces or commas, for example:

h->Draw("E1 SAME");

The options are not case sensitive:

h->Draw("e1 same");

The default drawing option can be set with TH1::SetOption and retrieve using TH1::GetOption:

root [0] h->Draw();          // Draw "h" using the standard histogram representation.
root [1] h->Draw("E");       // Draw "h" using error bars
root [3] h->SetOption("E");  // Change the default drawing option for "h"
root [4] h->Draw();          // Draw "h" using error bars
root [5] h->GetOption();     // Retrieve the default drawing option for "h"
(const Option_t* 0xa3ff948)"E"

Options supported for 1D and 2D histograms

Option Description
"E" Draw error bars.
"AXIS" Draw only axis.
"AXIG" Draw only grid (if the grid is requested).
"HIST" When an histogram has errors it is visualized by default with error bars. To visualize it without errors use the option "HIST" together with the required option (eg "hist same c"). The "HIST" option can also be used to plot only the histogram and not the associated function(s).
"FUNC" When an histogram has a fitted function, this option allows to draw the fit result only.
"SAME" Superimpose on previous picture in the same pad.
"SAMES" Same as "SAME" and draw the statistics box
"PFC" Palette Fill Color: histogram's fill color is taken in the current palette.
"PLC" Palette Line Color: histogram's line color is taken in the current palette.
"PMC" Palette Marker Color: histogram's marker color is taken in the current palette.
"LEGO" Draw a lego plot with hidden line removal.
"LEGO1" Draw a lego plot with hidden surface removal.
"LEGO2" Draw a lego plot using colors to show the cell contents When the option "0" is used with any LEGO option, the empty bins are not drawn.
"LEGO3" Draw a lego plot with hidden surface removal, like LEGO1 but the border lines of each lego-bar are not drawn.
"LEGO4" Draw a lego plot with hidden surface removal, like LEGO1 but without the shadow effect on each lego-bar.
"TEXT" Draw bin contents as text (format set via gStyle->SetPaintTextFormat).
"TEXTnn" Draw bin contents as text at angle nn (0 < nn < 90).
"X+" The X-axis is drawn on the top side of the plot.
"Y+" The Y-axis is drawn on the right side of the plot.
"MIN0" Set minimum value for the Y axis to 0, equivalent to gStyle->SetHistMinimumZero().

Options supported for 1D histograms

Option Description
" " Default.
"AH" Draw histogram without axis. "A" can be combined with any drawing option. For instance, "AC" draws the histogram as a smooth Curve without axis.
"][" When this option is selected the first and last vertical lines of the histogram are not drawn.
"B" Bar chart option.
"BAR" Like option "B", but bars can be drawn with a 3D effect.
"HBAR" Like option "BAR", but bars are drawn horizontally.
"C" Draw a smooth Curve through the histogram bins.
"E0" Draw error bars. Markers are drawn for bins with 0 contents.
"E1" Draw error bars with perpendicular lines at the edges.
"E2" Draw error bars with rectangles.
"E3" Draw a fill area through the end points of the vertical error bars.
"E4" Draw a smoothed filled area through the end points of the error bars.
"E5" Like E3 but ignore the bins with 0 contents.
"E6" Like E4 but ignore the bins with 0 contents.
"X0" When used with one of the "E" option, it suppress the error bar along X as gStyle->SetErrorX(0) would do.
"L" Draw a line through the bin contents.
"P" Draw current marker at each bin except empty bins.
"P0" Draw current marker at each bin including empty bins.
"PIE" Draw histogram as a Pie Chart.
"*H" Draw histogram with a * at each bin.
"LF2" Draw histogram like with option "L" but with a fill area. Note that "L" draws also a fill area if the hist fill color is set but the fill area corresponds to the histogram contour.

Options supported for 2D histograms

Option Description
" " Default (scatter plot).
"ARR" Arrow mode. Shows gradient between adjacent cells.
"BOX" A box is drawn for each cell with surface proportional to the content's absolute value. A negative content is marked with a X.
"BOX1" A button is drawn for each cell with surface proportional to content's absolute value. A sunken button is drawn for negative values a raised one for positive.
"COL" A box is drawn for each cell with a color scale varying with contents. All the none empty bins are painted. Empty bins are not painted unless some bins have a negative content because in that case the null bins might be not empty. TProfile2D histograms are handled differently because, for this type of 2D histograms, it is possible to know if an empty bin has been filled or not. So even if all the bins' contents are positive some empty bins might be painted. And vice versa, if some bins have a negative content some empty bins might be not painted.
"COLZ" Same as "COL". In addition the color palette is also drawn.
"COL2" Alternative rendering algorithm to "COL". Can significantly improve rendering performance for large, non-sparse 2-D histograms.
"COLZ2" Same as "COL2". In addition the color palette is also drawn.
"Z CJUST" In combination with colored options "COL","CONT0" etc: Justify labels in the color palette at color boudaries. For more details see TPaletteAxis
"CANDLE" Draw a candle plot along X axis.
"CANDLEX" Same as "CANDLE".
"CANDLEY" Draw a candle plot along Y axis.
"CANDLEXn" Draw a candle plot along X axis. Different candle-styles with n from 1 to 6.
"CANDLEYn" Draw a candle plot along Y axis. Different candle-styles with n from 1 to 6.
"VIOLIN" Draw a violin plot along X axis.
"VIOLINX" Same as "VIOLIN".
"VIOLINY" Draw a violin plot along Y axis.
"VIOLINXn" Draw a violin plot along X axis. Different violin-styles with n being 1 or 2.
"VIOLINYn" Draw a violin plot along Y axis. Different violin-styles with n being 1 or 2.
"CONT" Draw a contour plot (same as CONT0).
"CONT0" Draw a contour plot using surface colors to distinguish contours.
"CONT1" Draw a contour plot using line styles to distinguish contours.
"CONT2" Draw a contour plot using the same line style for all contours.
"CONT3" Draw a contour plot using fill area colors.
"CONT4" Draw a contour plot using surface colors (SURF option at theta = 0).
"CONT5" (TGraph2D only) Draw a contour plot using Delaunay triangles.
"LIST" Generate a list of TGraph objects for each contour.
"CYL" Use Cylindrical coordinates. The X coordinate is mapped on the angle and the Y coordinate on the cylinder length.
"POL" Use Polar coordinates. The X coordinate is mapped on the angle and the Y coordinate on the radius.
"SAME0" Same as "SAME" but do not use the z-axis range of the first plot.
"SAMES0" Same as "SAMES" but do not use the z-axis range of the first plot.
"SPH" Use Spherical coordinates. The X coordinate is mapped on the latitude and the Y coordinate on the longitude.
"PSR" Use PseudoRapidity/Phi coordinates. The X coordinate is mapped on Phi.
"SURF" Draw a surface plot with hidden line removal.
"SURF1" Draw a surface plot with hidden surface removal.
"SURF2" Draw a surface plot using colors to show the cell contents.
"SURF3" Same as SURF with in addition a contour view drawn on the top.
"SURF4" Draw a surface using Gouraud shading.
"SURF5" Same as SURF3 but only the colored contour is drawn. Used with option CYL, SPH or PSR it allows to draw colored contours on a sphere, a cylinder or a in pseudo rapidity space. In cartesian or polar coordinates, option SURF3 is used.
"LEGO9" Draw the 3D axis only. Mainly needed for internal use
"FB" With LEGO or SURFACE, suppress the Front-Box.
"BB" With LEGO or SURFACE, suppress the Back-Box.
"A" With LEGO or SURFACE, suppress the axis.
"SCAT" Draw a scatter-plot (default).
"[cutg]" Draw only the sub-range selected by the TCutG named "cutg".

Options supported for 3D histograms

Option Description
" " Default (scatter plot).
"ISO" Draw a Gouraud shaded 3d iso surface through a 3d histogram. It paints one surface at the value computed as follow: SumOfWeights/(NbinsX*NbinsY*NbinsZ).
"BOX" Draw a for each cell with volume proportional to the content's absolute value. An hidden line removal algorithm is used
"BOX1" Same as BOX but an hidden surface removal algorithm is used
"BOX2" The boxes' colors are picked in the current palette according to the bins' contents
"BOX2Z" Same as "BOX2". In addition the color palette is also drawn.
"BOX3" Same as BOX1, but the border lines of each lego-bar are not drawn.
"LEGO" Same as BOX.

Options supported for histograms' stacks (THStack)

Option Description
" " Default, the histograms are drawn on top of each other (as lego plots for 2D histograms).
"NOSTACK" Histograms in the stack are all paint in the same pad as if the option SAME had been specified.
"NOSTACKB" Histograms are drawn next to each other as bar charts.
"PADS" The current pad/canvas is subdivided into a number of pads equal to the number of histograms in the stack and each histogram is paint into a separate pad.
"PFC" Palette Fill Color: stack's fill color is taken in the current palette.
"PLC" Palette Line Color: stack's line color is taken in the current palette.
"PMC" Palette Marker Color: stack's marker color is taken in the current palette.

Setting the Style

Histograms use the current style (gStyle). When one changes the current style and would like to propagate the changes to the histogram, TH1::UseCurrentStyle should be called. Call UseCurrentStyle on each histogram is needed.

To force all the histogram to use the current style use:

gROOT->ForceStyle();

All the histograms read after this call will use the current style.

Setting line, fill, marker, and text attributes

The histogram classes inherit from the attribute classes: TAttLine, TAttFill and TAttMarker. See the description of these classes for the list of options.

Setting Tick marks on the histogram axis

The TPad::SetTicks method specifies the type of tick marks on the axis. If tx = gPad->GetTickx() and ty = gPad->GetTicky() then:

tx = 1;   tick marks on top side are drawn (inside)
tx = 2;   tick marks and labels on top side are drawn
ty = 1;   tick marks on right side are drawn (inside)
ty = 2;   tick marks and labels on right side are drawn

By default only the left Y axis and X bottom axis are drawn (tx = ty = 0)

TPad::SetTicks(tx,ty) allows to set these options. See also The TAxis functions to set specific axis attributes.

In case multiple color filled histograms are drawn on the same pad, the fill area may hide the axis tick marks. One can force a redraw of the axis over all the histograms by calling:

gPad->RedrawAxis();

Giving titles to the X, Y and Z axis

h->GetXaxis()->SetTitle("X axis title");
h->GetYaxis()->SetTitle("Y axis title");

The histogram title and the axis titles can be any TLatex string. The titles are part of the persistent histogram.

The option "SAME"

By default, when an histogram is drawn, the current pad is cleared before drawing. In order to keep the previous drawing and draw on top of it the option SAME should be use. The histogram drawn with the option SAME uses the coordinates system available in the current pad.

This option can be used alone or combined with any valid drawing option but some combinations must be use with care.

Limitations

  • It does not work when combined with the LEGO and SURF options unless the histogram plotted with the option SAME has exactly the same ranges on the X, Y and Z axis as the currently drawn histogram. To superimpose lego plots histograms' stacks should be used.

Colors automatically picked in palette

Since
ROOT version 6.09/01

When several histograms are painted in the same canvas thanks to the option "SAME" or via a THStack it might be useful to have an easy and automatic way to choose their color. The simplest way is to pick colors in the current active color palette. Palette coloring for histogram is activated thanks to the options PFC (Palette Fill Color), PLC (Palette Line Color) and PMC (Palette Marker Color). When one of these options is given to TH1::Draw the histogram get its color from the current color palette defined by gStyle->SetPalette(…). The color is determined according to the number of objects having palette coloring in the current pad.

void histpalettecolor()
{
auto C = new TCanvas();
auto h1 = new TH1F ("h1","Histogram drawn with full circles",100,-4,4);
auto h2 = new TH1F ("h2","Histogram drawn with full squares",100,-4,4);
auto h3 = new TH1F ("h3","Histogram drawn with full triangles up",100,-4,4);
auto h4 = new TH1F ("h4","Histogram drawn with full triangles down",100,-4,4);
auto h5 = new TH1F ("h5","Histogram drawn with empty circles",100,-4,4);
TRandom3 rng;
Double_t px,py;
for (Int_t i = 0; i < 25000; i++) {
rng.Rannor(px,py);
h1->Fill(px,10.);
h2->Fill(px, 8.);
h3->Fill(px, 6.);
h4->Fill(px, 4.);
h5->Fill(px, 2.);
}
h2->SetMarkerStyle(kFullSquare);
h3->SetMarkerStyle(kFullTriangleUp);
h4->SetMarkerStyle(kFullTriangleDown);
h5->SetMarkerStyle(kOpenCircle);
h1->Draw("PLC PMC");
h2->Draw("SAME PLC PMC");
h3->Draw("SAME PLC PMC");
h4->Draw("SAME PLC PMC");
h5->Draw("SAME PLC PMC");
gPad->BuildLegend();
}
int Int_t
Definition: RtypesCore.h:41
const Bool_t kFALSE
Definition: RtypesCore.h:88
double Double_t
Definition: RtypesCore.h:55
@ kFullTriangleDown
Definition: TAttMarker.h:49
@ kFullSquare
Definition: TAttMarker.h:48
@ kFullTriangleUp
Definition: TAttMarker.h:48
@ kFullCircle
Definition: TAttMarker.h:48
@ kOpenCircle
Definition: TAttMarker.h:49
R__EXTERN TStyle * gStyle
Definition: TStyle.h:407
#define gPad
Definition: TVirtualPad.h:286
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
Definition: TAttMarker.h:40
The Canvas class.
Definition: TCanvas.h:31
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:571
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition: TH1.cxx:3275
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2998
Random number generator class based on M.
Definition: TRandom3.h:27
virtual void Rannor(Float_t &a, Float_t &b)
Return 2 numbers distributed following a gaussian with mean=0 and sigma=1.
Definition: TRandom.cxx:489
void SetOptTitle(Int_t tit=1)
Definition: TStyle.h:313
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:1450
TH1F * h1
Definition: legend1.C:5
static double C[]
void thstackpalettecolor()
{
auto hs = new THStack("hs","Stacked 1D histograms colored using kOcean palette");
// Create three 1-d histograms and add them in the stack
auto h1st = new TH1F("h1st","test hstack",100,-4,4);
h1st->FillRandom("gaus",20000);
hs->Add(h1st);
auto h2st = new TH1F("h2st","test hstack",100,-4,4);
h2st->FillRandom("gaus",15000);
hs->Add(h2st);
auto h3st = new TH1F("h3st","test hstack",100,-4,4);
h3st->FillRandom("gaus",10000);
hs->Add(h3st);
// draw the stack
hs->Draw("pfc nostack");
}
@ kOcean
Definition: TColor.h:110
The Histogram stack class.
Definition: THStack.h:31
void SetPalette(Int_t ncolors=kBird, Int_t *colors=0, Float_t alpha=1.)
See TColor::SetPalette.
Definition: TStyle.cxx:1643
void thstack2palettecolor () {
auto h1 = new TH2F("h1","h1",20,0,6,20,-4,4);
auto h2 = new TH2F("h2","h1",20,0,6,20,-4,4);
auto h3 = new TH2F("h3","h1",20,0,6,20,-4,4);
auto h4 = new TH2F("h4","h1",20,0,6,20,-4,4);
auto h5 = new TH2F("h5","h1",20,0,6,20,-4,4);
h2->Fill(2.,0.,5);
h3->Fill(3.,0.,10);
h4->Fill(4.,0.,15);
h5->Fill(5.,0.,20);
auto hs = new THStack("hs","Test of palette colored lego stack");
hs->Add(h1);
hs->Add(h2);
hs->Add(h3);
hs->Add(h4);
hs->Add(h5);
hs->Draw("0lego1 PFC");
}
2-D histogram with a float per channel (see TH1 documentation)}
Definition: TH2.h:251

Superimposing two histograms with different scales in the same pad

The following example creates two histograms, the second histogram is the bins integral of the first one. It shows a procedure to draw the two histograms in the same pad and it draws the scale of the second histogram using a new vertical axis on the right side. See also the tutorial transpad.C for a variant of this example.

{
auto c1 = new TCanvas("c1","c1",600,400);
// create/fill draw h1
auto h1 = new TH1F("h1","Superimposing two histograms with different scales",100,-3,3);
Int_t i;
for (i=0;i<10000;i++) h1->Fill(gRandom->Gaus(0,1));
h1->Draw();
c1->Update();
// create hint1 filled with the bins integral of h1
auto hint1 = new TH1F("hint1","h1 bins integral",100,-3,3);
float sum = 0.f;
for (i=1;i<=100;i++) {
hint1->SetBinContent(i,sum);
}
// scale hint1 to the pad coordinates
float rightmax = 1.1*hint1->GetMaximum();
float scale = gPad->GetUymax()/rightmax;
hint1->SetLineColor(kRed);
hint1->Scale(scale);
hint1->Draw("same");
// draw an axis on the right side
auto axis = new TGaxis(gPad->GetUxmax(),gPad->GetUymin(),
gPad->GetUxmax(), gPad->GetUymax(),0,rightmax,510,"+L");
axis->SetLineColor(kRed);
axis->SetTextColor(kRed);
axis->Draw();
}
@ kRed
Definition: Rtypes.h:64
R__EXTERN TRandom * gRandom
Definition: TRandom.h:62
The axis painter class.
Definition: TGaxis.h:24
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition: TH1.cxx:4899
virtual Double_t Gaus(Double_t mean=0, Double_t sigma=1)
Samples a random number from the standard Normal (Gaussian) Distribution with the given mean and sigm...
Definition: TRandom.cxx:263
return c1
Definition: legend1.C:41
static long int sum(long int i)
Definition: Factory.cxx:2276

Statistics Display

The type of information shown in the histogram statistics box can be selected with:

gStyle->SetOptStat(mode);

The mode has up to nine digits that can be set to on (1 or 2), off (0).

mode = ksiourmen  (default = 000001111)
k = 1;  kurtosis printed
k = 2;  kurtosis and kurtosis error printed
s = 1;  skewness printed
s = 2;  skewness and skewness error printed
i = 1;  integral of bins printed
i = 2;  integral of bins with option "width" printed
o = 1;  number of overflows printed
u = 1;  number of underflows printed
r = 1;  standard deviation printed
r = 2;  standard deviation and standard deviation error printed
m = 1;  mean value printed
m = 2;  mean and mean error values printed
e = 1;  number of entries printed
n = 1;  name of histogram is printed

For example:

gStyle->SetOptStat(11);

displays only the name of histogram and the number of entries, whereas:

gStyle->SetOptStat(1101);

displays the name of histogram, mean value and standard deviation.

WARNING 1: never do:

gStyle->SetOptStat(0001111);

but instead do:

gStyle->SetOptStat(1111);

because 0001111 will be taken as an octal number!

WARNING 2: for backward compatibility with older versions

gStyle->SetOptStat(1);

is taken as:

gStyle->SetOptStat(1111)

To print only the name of the histogram do:

gStyle->SetOptStat(1000000001);

NOTE that in case of 2D histograms, when selecting only underflow (10000) or overflow (100000), the statistics box will show all combinations of underflow/overflows and not just one single number.

The parameter mode can be any combination of the letters kKsSiIourRmMen

k :  kurtosis printed
K :  kurtosis and kurtosis error printed
s :  skewness printed
S :  skewness and skewness error printed
i :  integral of bins printed
I :  integral of bins with option "width" printed
o :  number of overflows printed
u :  number of underflows printed
r :  standard deviation printed
R :  standard deviation and standard deviation error printed
m :  mean value printed
M :  mean value mean error values printed
e :  number of entries printed
n :  name of histogram is printed

For example, to print only name of histogram and number of entries do:

gStyle->SetOptStat("ne");

To print only the name of the histogram do:

gStyle->SetOptStat("n");

The default value is:

gStyle->SetOptStat("nemr");

When a histogram is painted, a TPaveStats object is created and added to the list of functions of the histogram. If a TPaveStats object already exists in the histogram list of functions, the existing object is just updated with the current histogram parameters.

Once a histogram is painted, the statistics box can be accessed using h->FindObject("stats"). In the command line it is enough to do:

Root > h->Draw()
Root > TPaveStats *st = (TPaveStats*)h->FindObject("stats")

because after h->Draw() the histogram is automatically painted. But in a script file the painting should be forced using gPad->Update() in order to make sure the statistics box is created:

h->Draw();
gPad->Update();
TPaveStats *st = (TPaveStats*)h->FindObject("stats");

Without gPad->Update() the line h->FindObject("stats") returns a null pointer.

When a histogram is drawn with the option SAME, the statistics box is not drawn. To force the statistics box drawing with the option SAME, the option SAMES must be used. If the new statistics box hides the previous statistics box, one can change its position with these lines (h being the pointer to the histogram):

Root > TPaveStats *st = (TPaveStats*)h->FindObject("stats")
Root > st->SetX1NDC(newx1); //new x start position
Root > st->SetX2NDC(newx2); //new x end position

To change the type of information for an histogram with an existing TPaveStats one should do:

st->SetOptStat(mode);

Where mode has the same meaning than when calling gStyle->SetOptStat(mode) (see above).

One can delete the statistics box for a histogram TH1* h with:

h->SetStats(0)

and activate it again with:

h->SetStats(1).

Labels used in the statistics box ("Mean", "Std Dev", ...) can be changed from $ROOTSYS/etc/system.rootrc or .rootrc (look for the string Hist.Stats.).

Fit Statistics

The type of information about fit parameters printed in the histogram statistics box can be selected via the parameter mode. The parameter mode can be = pcev (default = 0111)

p = 1;  print Probability
c = 1;  print Chisquare/Number of degrees of freedom
e = 1;  print errors (if e=1, v must be 1)
v = 1;  print name/values of parameters

Example:

gStyle->SetOptFit(1011);

print fit probability, parameter names/values and errors.

  1. When v" = 1</tt> is specified, only the non-fixed parameters are shown. 2. When <tt>v" = 2 all parameters are shown.

Note: gStyle->SetOptFit(1) means "default value", so it is equivalent to gStyle->SetOptFit(111)

The error bars options

Option Description
"E" Default. Shows only the error bars, not a marker.
"E1" Small lines are drawn at the end of the error bars.
"E2" Error rectangles are drawn.
"E3" A filled area is drawn through the end points of the vertical error bars.
"E4" A smoothed filled area is drawn through the end points of the vertical error bars.
"E0" Draw also bins with null contents.
{
auto c1 = new TCanvas("c1","c1",600,400);
auto he = new TH1F("he","Distribution drawn with error bars (option E1) ",100,-3,3);
for (int i=0; i<10000; i++) he->Fill(gRandom->Gaus(0,1));
he->SetMarkerStyle(20);
he->Draw("E1");
}
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:1152
void SetErrorX(Float_t errorx=0.5)
Definition: TStyle.h:319

The options "E3" and "E4" draw an error band through the end points of the vertical error bars. With "E4" the error band is smoothed. Because of the smoothing algorithm used some artefacts may appear at the end of the band like in the following example. In such cases "E3" should be used instead of "E4".

{
auto ce4 = new TCanvas("ce4","ce4",600,400);
ce4->Divide(2,1);
auto he4 = new TH1F("he4","Distribution drawn with option E4",100,-3,3);
Int_t i;
for (i=0;i<10000;i++) he4->Fill(gRandom->Gaus(0,1));
he4->SetFillColor(kRed);
he4->GetXaxis()->SetRange(40,48);
ce4->cd(1);
he4->Draw("E4");
ce4->cd(2);
auto he3 = (TH1F*)he4->DrawClone("E3");
he3->SetTitle("Distribution drawn option E3");
}

2D histograms can be drawn with error bars as shown is the following example:

{
auto c2e = new TCanvas("c2e","c2e",600,400);
auto h2e = new TH2F("h2e","TH2 drawn with option E",40,-4,4,40,-20,20);
float px, py;
for (Int_t i = 0; i < 25000; i++) {
gRandom->Rannor(px,py);
h2e->Fill(px,5*py);
}
h2e->Draw("E");
}

The bar chart option

The option "B" allows to draw simple vertical bar charts. The bar width is controlled with TH1::SetBarWidth(), and the bar offset within the bin, with TH1::SetBarOffset(). These two settings are useful to draw several histograms on the same plot as shown in the following example:

{
int i;
const Int_t nx = 8;
string os_X[nx] = {"8","32","128","512","2048","8192","32768","131072"};
float d_35_0[nx] = {0.75, -3.30, -0.92, 0.10, 0.08, -1.69, -1.29, -2.37};
float d_35_1[nx] = {1.01, -3.02, -0.65, 0.37, 0.34, -1.42, -1.02, -2.10};
auto cb = new TCanvas("cb","cb",600,400);
cb->SetGrid();
auto h1b = new TH1F("h1b","Option B example",nx,0,nx);
h1b->SetFillColor(4);
h1b->SetBarWidth(0.4);
h1b->SetBarOffset(0.1);
h1b->SetStats(0);
h1b->SetMinimum(-5);
h1b->SetMaximum(5);
for (i=1; i<=nx; i++) {
h1b->SetBinContent(i, d_35_0[i-1]);
h1b->GetXaxis()->SetBinLabel(i,os_X[i-1].c_str());
}
h1b->Draw("b");
auto h2b = new TH1F("h2b","h2b",nx,0,nx);
h2b->SetFillColor(38);
h2b->SetBarWidth(0.4);
h2b->SetBarOffset(0.5);
h2b->SetStats(0);
for (i=1;i<=nx;i++) h2b->SetBinContent(i, d_35_1[i-1]);
h2b->Draw("b same");
}
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:1100

The "BAR" and "HBAR" options

When the option bar or hbar is specified, a bar chart is drawn. A vertical bar-chart is drawn with the options bar, bar0, bar1, bar2, bar3, bar4. An horizontal bar-chart is drawn with the options hbar, hbar0, hbar1, hbar2, hbar3, hbar4 (hbars.C).

  • The bar is filled with the histogram fill color.
  • The left side of the bar is drawn with a light fill color.
  • The right side of the bar is drawn with a dark fill color.
  • The percentage of the bar drawn with either the light or dark color is:
    • 0% for option "(h)bar" or "(h)bar0"
    • 10% for option "(h)bar1"
    • 20% for option "(h)bar2"
    • 30% for option "(h)bar3"
    • 40% for option "(h)bar4"

When an histogram has errors the option "HIST" together with the (h)bar option.

TCanvas *hbars() {
// Try to open first the file cernstaff.root in tutorials/tree directory
TString filedir = gROOT->GetTutorialDir();
filedir += TString("/tree/");
TString filename = "cernstaff.root";
bool fileNotFound = gSystem->AccessPathName(filename); // note opposite return code
// If file is not found try to generate it uing the macro tree/cernbuild.C
if (fileNotFound) {
TString macroName = filedir + "cernbuild.C";
if (!gInterpreter->IsLoaded(macroName)) gInterpreter->LoadMacro(macroName);
gROOT->ProcessLineFast("cernbuild()");
}
TFile * f = TFile::Open(filename);
if (!f) {
Error("hbars","file cernstaff.root not found");
return 0;
}
TTree *T = (TTree*)f->Get("T");
if (!T) {
Error("hbars","Tree T is not present in file %s",f->GetName() );
return 0;
}
T->SetFillColor(45);
TCanvas *c1 = new TCanvas("c1","histograms with bars",700,800);
c1->SetFillColor(42);
c1->Divide(1,2);
// Horizontal bar chart
c1->cd(1); gPad->SetGrid(); gPad->SetLogx(); gPad->SetFrameFillColor(33);
T->Draw("Nation","","hbar2");
// Vertical bar chart
c1->cd(2); gPad->SetGrid(); gPad->SetFrameFillColor(33);
T->Draw("Division>>hDiv","","goff");
TH1F *hDiv = (TH1F*)gDirectory->Get("hDiv");
hDiv->SetStats(0);
TH1F *hDivFR = (TH1F*)hDiv->Clone("hDivFR");
T->Draw("Division>>hDivFR","Nation==\"FR\"","goff");
hDiv->SetBarWidth(0.45);
hDiv->SetBarOffset(0.1);
hDiv->SetFillColor(49);
TH1 *h1 = hDiv->DrawCopy("bar2");
hDivFR->SetBarWidth(0.4);
hDivFR->SetBarOffset(0.55);
hDivFR->SetFillColor(50);
TH1 *h2 = hDivFR->DrawCopy("bar2,same");
TLegend *legend = new TLegend(0.55,0.65,0.76,0.82);
legend->AddEntry(h1,"All nations","f");
legend->AddEntry(h2,"French only","f");
legend->Draw();
c1->cd();
delete f;
return c1;
}
#define f(i)
Definition: RSha256.hxx:104
#define gDirectory
Definition: TDirectory.h:223
#define gInterpreter
Definition: TInterpreter.h:555
#define gROOT
Definition: TROOT.h:415
R__EXTERN TSystem * gSystem
Definition: TSystem.h:560
virtual void SetFillColor(Color_t fcolor)
Set the fill area color.
Definition: TAttFill.h:37
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format.
Definition: TFile.h:48
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=ROOT::RCompressionSetting::EDefaults::kUseCompiledDefault, Int_t netopt=0)
Create / open a file.
Definition: TFile.cxx:3923
The TH1 histogram class.
Definition: TH1.h:56
TObject * Clone(const char *newname=0) const
Make a complete copy of the underlying object.
Definition: TH1.cxx:2665
virtual void SetBarWidth(Float_t width=0.5)
Definition: TH1.h:356
virtual TH1 * DrawCopy(Option_t *option="", const char *name_postfix="_copy") const
Copy this histogram and Draw in the current pad.
Definition: TH1.cxx:3045
virtual void SetBarOffset(Float_t offset=0.25)
Definition: TH1.h:355
virtual void SetStats(Bool_t stats=kTRUE)
Set statistics option on/off.
Definition: TH1.cxx:8434
This class displays a legend box (TPaveText) containing several legend entries.
Definition: TLegend.h:23
TLegendEntry * AddEntry(const TObject *obj, const char *label="", Option_t *option="lpf")
Add a new entry to this legend.
Definition: TLegend.cxx:330
virtual void Draw(Option_t *option="")
Draw this legend with its current attributes.
Definition: TLegend.cxx:423
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:880
Basic string class.
Definition: TString.h:131
virtual Bool_t AccessPathName(const char *path, EAccessMode mode=kFileExists)
Returns FALSE if one can access a file using the specified access mode.
Definition: TSystem.cxx:1287
A TTree represents a columnar dataset.
Definition: TTree.h:72
double T(double x)
Definition: ChebyshevPol.h:34

To control the bar width (default is the bin width) TH1::SetBarWidth() should be used.

To control the bar offset (default is 0) TH1::SetBarOffset() should be used.

These two parameters are useful when several histograms are plotted using the option SAME. They allow to plot the histograms next to each other.

The SCATter plot option (default for 2D histograms)

For each cell (i,j) a number of points proportional to the cell content is drawn. A maximum of kNMAX points per cell is drawn. If the maximum is above kNMAX contents are normalized to kNMAX (kNMAX=2000). If option is of the form scat=ff, (eg scat=1.8, scat=1e-3), then ff is used as a scale factor to compute the number of dots. scat=1 is the default.

By default the scatter plot is painted with a "dot marker" which not scalable (see the TAttMarker documentation). To change the marker size, a scalable marker type should be used. For instance a circle (marker style 20).

{
auto c1 = new TCanvas("c1","c1",600,400);
auto hscat = new TH2F("hscat","Option SCATter example (default for 2D histograms) ",40,-4,4,40,-20,20);
float px, py;
for (Int_t i = 0; i < 25000; i++) {
gRandom->Rannor(px,py);
hscat->Fill(px,5*py);
hscat->Fill(3+0.5*px,2*py-10.);
}
hscat->Draw("scat=0.5");
}

The ARRow option

Shows gradient between adjacent cells. For each cell (i,j) an arrow is drawn The orientation of the arrow follows the cell gradient.

{
auto c1 = new TCanvas("c1","c1",600,400);
auto harr = new TH2F("harr","Option ARRow example",20,-4,4,20,-20,20);
harr->SetLineColor(kRed);
float px, py;
for (Int_t i = 0; i < 25000; i++) {
gRandom->Rannor(px,py);
harr->Fill(px,5*py);
harr->Fill(3+0.5*px,2*py-10.,0.1);
}
harr->Draw("ARR");
}
Since
ROOT version 6.17/01

The option ARR can be combined with the option COL or COLZ.

{
auto c1 = new TCanvas("c1","c1",600,400);
auto harr = new TH2F("harr","Option ARR + COLZ example",20,-4,4,20,-20,20);
harr->SetStats(0);
float px, py;
for (Int_t i = 0; i < 25000; i++) {
gRandom->Rannor(px,py);
harr->Fill(px,5*py);
harr->Fill(3+0.5*px,2*py-10.,0.1);
}
harr->Draw("ARR COLZ");
}

The BOX option

For each cell (i,j) a box is drawn. The size (surface) of the box is proportional to the absolute value of the cell content. The cells with a negative content are drawn with a X on top of the box.

{
auto c1 = new TCanvas("c1","c1",600,400);
auto hbox = new TH2F("hbox","Option BOX example",3,0,3,3,0,3);
hbox->SetFillColor(42);
hbox->Fill(0.5, 0.5, 1.);
hbox->Fill(0.5, 1.5, 4.);
hbox->Fill(0.5, 2.5, 3.);
hbox->Fill(1.5, 0.5, 2.);
hbox->Fill(1.5, 1.5, 12.);
hbox->Fill(1.5, 2.5, -6.);
hbox->Fill(2.5, 0.5, -4.);
hbox->Fill(2.5, 1.5, 6.);
hbox->Fill(2.5, 2.5, 0.5);
hbox->Draw("BOX");
}

With option BOX1 a button is drawn for each cell with surface proportional to content's absolute value. A sunken button is drawn for negative values a raised one for positive.

{
auto c1 = new TCanvas("c1","c1",600,400);
auto hbox1 = new TH2F("hbox1","Option BOX1 example",3,0,3,3,0,3);
hbox1->SetFillColor(42);
hbox1->Fill(0.5, 0.5, 1.);
hbox1->Fill(0.5, 1.5, 4.);
hbox1->Fill(0.5, 2.5, 3.);
hbox1->Fill(1.5, 0.5, 2.);
hbox1->Fill(1.5, 1.5, 12.);
hbox1->Fill(1.5, 2.5, -6.);
hbox1->Fill(2.5, 0.5, -4.);
hbox1->Fill(2.5, 1.5, 6.);
hbox1->Fill(2.5, 2.5, 0.5);
hbox1->Draw("BOX1");
}

When the option SAME (or "SAMES") is used with the option BOX, the boxes' sizes are computing taking the previous plots into account. The range along the Z axis is imposed by the first plot (the one without option SAME); therefore the order in which the plots are done is relevant.

{
auto c1 = new TCanvas("c1","c1",600,400);
auto hb1 = new TH2F("hb1","Example of BOX plots with option SAME ",40,-3,3,40,-3,3);
auto hb2 = new TH2F("hb2","hb2",40,-3,3,40,-3,3);
auto hb3 = new TH2F("hb3","hb3",40,-3,3,40,-3,3);
auto hb4 = new TH2F("hb4","hb4",40,-3,3,40,-3,3);
for (Int_t i=0;i<1000;i++) {
double x,y;
if (x>0 && y>0) hb1->Fill(x,y,4);
if (x<0 && y<0) hb2->Fill(x,y,3);
if (x>0 && y<0) hb3->Fill(x,y,2);
if (x<0 && y>0) hb4->Fill(x,y,1);
}
hb1->SetFillColor(1);
hb2->SetFillColor(2);
hb3->SetFillColor(3);
hb4->SetFillColor(4);
hb1->Draw("box");
hb2->Draw("box same");
hb3->Draw("box same");
hb4->Draw("box same");
}
Double_t y[n]
Definition: legend1.C:17
Double_t x[n]
Definition: legend1.C:17
Since
ROOT version 6.17/01:

Sometimes the change of the range of the Z axis is unwanted, in which case, one can use SAME0 (or SAMES0) option to opt out of this change.

{
auto h2 = new TH2F("h2"," ",10,0,10,10,20,30);
auto hf = (TH2F*)h2->Clone("hf");
hf->SetBit(TH1::kNoStats);
h2->Fill(5,22);
h2->Fill(5,23);
h2->Fill(6,22);
h2->Fill(6,23);
hf->Fill(6,23);
hf->Fill(6,23);
hf->Fill(6,23);
hf->Fill(6,23);
hf->Fill(5,23);
auto hf_copy1 = hf->Clone("hf_copy1");
auto lt = new TLatex();
auto cx = new TCanvas(); cx->Divide(2,1);
cx->cd(1);
h2->Draw("box");
hf->Draw("text colz same");
lt->DrawLatexNDC(0.3,0.5,"SAME");
cx->cd(2);
h2->Draw("box");
hf_copy1->Draw("text colz same0");
lt->DrawLatexNDC(0.3,0.5,"SAME0");
}
@ kNoStats
don't draw stats box
Definition: TH1.h:160
To draw Mathematical Formula.
Definition: TLatex.h:18
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition: TObject.cxx:694

The COLor option

For each cell (i,j) a box is drawn with a color proportional to the cell content.

The color table used is defined in the current style.

If the histogram's minimum and maximum are the same (flat histogram), the mapping on colors is not possible, therefore nothing is painted. To paint a flat histogram it is enough to set the histogram minimum (TH1::SetMinimum()) different from the bins' content.

The default number of color levels used to paint the cells is 20. It can be changed with TH1::SetContour() or TStyle::SetNumberContours(). The higher this number is, the smoother is the color change between cells.

The color palette in TStyle can be modified via gStyle->SetPalette().

All the non-empty bins are painted. Empty bins are not painted unless some bins have a negative content because in that case the null bins might be not empty.

TProfile2D histograms are handled differently because, for this type of 2D histograms, it is possible to know if an empty bin has been filled or not. So even if all the bins' contents are positive some empty bins might be painted. And vice versa, if some bins have a negative content some empty bins might be not painted.

Combined with the option COL, the option Z allows to display the color palette defined by gStyle->SetPalette().

In the following example, the histogram has only positive bins; the empty bins (containing 0) are not drawn.

{
auto c1 = new TCanvas("c1","c1",600,400);
auto hcol1 = new TH2F("hcol1","Option COLor example ",40,-4,4,40,-20,20);
float px, py;
for (Int_t i = 0; i < 25000; i++) {
gRandom->Rannor(px,py);
hcol1->Fill(px,5*py);
}
hcol1->Draw("COLZ");
}

In the first plot of following example, the histogram has some negative bins; the empty bins (containing 0) are drawn. In some cases one wants to not draw empty bins (containing 0) of histograms having a negative minimum. The option 1, used to produce the second plot in the following picture, allows to do that.

{
auto c1 = new TCanvas("c1","c1",600,600);
c1->Divide(1,2);
auto hcol23 = new TH2F("hcol23","Option COLZ example ",40,-4,4,40,-20,20);
auto hcol24 = new TH2F("hcol24","Option COLZ1 example ",40,-4,4,40,-20,20);
float px, py;
for (Int_t i = 0; i < 25000; i++) {
gRandom->Rannor(px,py);
hcol23->Fill(px,5*py);
hcol24->Fill(px,5*py);
}
hcol23->Fill(0.,0.,-200.);
hcol24->Fill(0.,0.,-200.);
c1->cd(1); hcol23->Draw("COLZ");
c1->cd(2); hcol24->Draw("COLZ1");
}

When the maximum of the histogram is set to a smaller value than the real maximum, the bins having a content between the new maximum and the real maximum are painted with the color corresponding to the new maximum.

When the minimum of the histogram is set to a greater value than the real minimum, the bins having a value between the real minimum and the new minimum are not drawn unless the option 0 is set.

The following example illustrates the option 0 combined with the option COL.

{
auto c1 = new TCanvas("c1","c1",600,600);
c1->Divide(1,2);
auto hcol21 = new TH2F("hcol21","Option COLZ",40,-4,4,40,-20,20);
auto hcol22 = new TH2F("hcol22","Option COLZ0",40,-4,4,40,-20,20);
float px, py;
for (Int_t i = 0; i < 25000; i++) {
gRandom->Rannor(px,py);
hcol21->Fill(px,5*py);
hcol22->Fill(px,5*py);
}
hcol21->SetBit(TH1::kNoStats);
hcol22->SetBit(TH1::kNoStats);
c1->cd(1); hcol21->Draw("COLZ");
c1->cd(2); hcol22->Draw("COLZ0");
hcol22->SetMaximum(100);
hcol22->SetMinimum(40);
}
Since
ROOT version 6.09/01:

When the option SAME (or "SAMES") is used with the option COL, the boxes' color are computing taking the previous plots into account. The range along the Z axis is imposed by the first plot (the one without option SAME); therefore the order in which the plots are done is relevant. Same as in the `BOX` option, one can use SAME0 (or SAMES0) to opt out of this imposition.

{
auto c = new TCanvas("c","Example of col plots with option SAME",200,10,700,500);
auto h1 = new TH2F("h1","h1",40,-3,3,40,-3,3);
auto h2 = new TH2F("h2","h2",40,-3,3,40,-3,3);
auto h3 = new TH2F("h3","h3",40,-3,3,40,-3,3);
auto h4 = new TH2F("h4","h4",40,-3,3,40,-3,3);
for (Int_t i=0;i<5000;i++) {
double x,y;
if(x>0 && y>0) h1->Fill(x,y,4);
if(x<0 && y<0) h2->Fill(x,y,3);
if(x>0 && y<0) h3->Fill(x,y,2);
if(x<0 && y>0) h4->Fill(x,y,1);
}
h1->Draw("colz");
h2->Draw("col same");
h3->Draw("col same");
h4->Draw("col same");
}
#define c(i)
Definition: RSha256.hxx:101

The option COL can be combined with the option POL:

{
auto c1 = new TCanvas("c1","c1",600,400);
auto hcol1 = new TH2F("hcol1","Option COLor combined with POL",40,-4,4,40,-4,4);
float px, py;
for (Int_t i = 0; i < 25000; i++) {
gRandom->Rannor(px,py);
hcol1->Fill(px,py);
}
hcol1->Draw("COLZPOL");
}
Since
ROOT version 6.07/03:

A second rendering technique is also available with the COL2 and COLZ2 options.

These options provide potential performance improvements compared to the standard COL option. The performance comparison of the COL2 to the COL option depends on the histogram and the size of the rendering region in the current pad. In general, a small (approx. less than 100 bins per axis), sparsely populated TH2 will render faster with the COL option.

However, for larger histograms (approx. more than 100 bins per axis) that are not sparse, the COL2 option will provide up to 20 times performance improvements. For example, a 1000x1000 bin TH2 that is not sparse will render an order of magnitude faster with the COL2 option.

The COL2 option will also scale its performance based on the size of the pixmap the histogram image is being rendered into. It also is much better optimized for sessions where the user is forwarding X11 windows through an ssh connection.

For the most part, the COL2 and COLZ2 options are a drop in replacement to the COL and COLZ options. There is one major difference and that concerns the treatment of bins with zero content. The COL2 and COLZ2 options color these bins the color of zero.

COL2 option renders the histogram as a bitmap. Therefore it cannot be saved in vector graphics file format like PostScript or PDF (an empty image will be generated). It can be saved only in bitmap files like PNG format for instance.

The CANDLE and VIOLIN options

The mechanism behind Candle plots and Violin plots is very similar. Because of this they are implemented in the same class TCandle. The keywords CANDLE or VIOLIN will initiate the drawing of the corresponding plots. Followed by the keyword the user can select a plot direction (X or V for vertical projections, or Y or H for horizontal projections) and/or predefined definitions (1-6 for candles, 1-2 for violins). The order doesn't matter. Default is X and 1.

Instead of using the predefined representations, the candle and violin parameters can be changed individually. In that case the option have the following form:

CANDLEX(<option-string>)
CANDLEY(<option-string>)
VIOLINX(<option-string>)
VIOLINY(<option-string>).

All zeros at the beginning of option-string can be omitted.

option-string consists eight values, defined as follow:

"CANDLEX(zhpawMmb)"

Where:

  • b = 0; no box drawn
  • b = 1; the box is drawn. As the candle-plot is also called a box-plot it makes sense in the very most cases to always draw the box
  • b = 2; draw a filled box with border
  • m = 0; no median drawn
  • m = 1; median is drawn as a line
  • m = 2; median is drawn with errors (notches)
  • m = 3; median is drawn as a circle
  • M = 0; no mean drawn
  • M = 1; mean is drawn as a dashed line
  • M = 3; mean is drawn as a circle
  • w = 0; no whisker drawn
  • w = 1; whisker is drawn to end of distribution.
  • w = 2; whisker is drawn to max 1.5*iqr
  • a = 0; no anchor drawn
  • a = 1; the anchors are drawn
  • p = 0; no points drawn
  • p = 1; only outliers are drawn
  • p = 2; all datapoints are drawn
  • p = 3: all datapoints are drawn scattered
  • h = 0; no histogram is drawn
  • h = 1; histogram at the left or bottom side is drawn
  • h = 2; histogram at the right or top side is drawn
  • h = 3; histogram at left and right or top and bottom (violin-style) is drawn
  • z = 0; no zero indicator line is drawn
  • z = 1; zero indicator line is drawn.

As one can see all individual options for both candle and violin plots can be accessed by this mechanism. In deed the keywords CANDLE(<option-string>) and VIOLIN(<option-string>) have the same meaning. So you can parametrise an option-string for a candle plot and use the keywords VIOLIN and vice versa, if you wish.

Using a logarithmic x- or y-axis is possible for candle and violin charts.

Since
ROOT version 6.11/01

a logarithmic z-axis is possible, too but will only affect violin charts of course.

The CANDLE option

A Candle plot (also known as a "box plot" or "whisker plot") was invented in 1977 by John Tukey. It is a convenient way to describe graphically a data distribution (D) with only five numbers:

  1. The minimum value of the distribution D (bottom or left whisker).
  2. The lower quartile (Q1): 25% of the data points in D are less than Q1 (bottom of the box).
  3. The median (M): 50% of the data points in D are less than M.
  4. The upper quartile (Q3): 75% of the data points in D are less than Q3 (top of the box).
  5. The maximum value of the distribution D (top or right whisker).

In this implementation a TH2 is considered as a collection of TH1 along X (option CANDLE or CANDLEX) or Y (option CANDLEY). Each TH1 is represented as one candle.

void candleplotwhiskers() {
auto c1 = new TCanvas("c1","Candle Presets",700,800);
c1->Divide(1,2);
auto rng = new TRandom();
auto h1 = new TH2I("h1","Gaus",100,-5,5,1,0,1);
auto h2 = new TH1I("h2","Gaus",100,-5,5);
h1->GetXaxis()->SetTitle("Standard deviation #sigma");
h2->GetXaxis()->SetTitle("Standard deviation #sigma");
h2->GetYaxis()->SetTitle("dN/d#sigma");
float myRand;
for (int i = 0; i < 100000; i++) {
myRand = rng->Gaus(0,1);
h1->Fill(myRand,0);
h2->Fill(myRand);
}
Double_t *q = new Double_t[3];
Double_t *p = new Double_t[3];
q[0] = 0.; q[1] = 0.; q[2] = 0.;
p[0] = 0.25; p[1] = 0.5; p[2] = 0.75;
h2->GetQuantiles(3,q,p);
cout << "Q1 (-25%): " << q[0] << " Median: " << q[1] << " Q3 (+25%): " << q[2] << endl;
double iqr = q[2]-q[0];
auto mygaus_1_middle = new TF1("mygaus_1_middle","gaus",q[0],q[2]);
auto mygaus_1_left = new TF1("mygaus_1_left","gaus",q[0]-1.5*iqr,q[0]);
mygaus_1_left->SetLineColor(kGreen);
auto mygaus_1_right = new TF1("mygaus_1_right","gaus",q[2],q[2]+1.5*iqr);
mygaus_1_right->SetLineColor(kGreen);
c1->cd(1);
h1->Draw("candley2 scat");
c1->cd(2);
h2->Draw("");
h2->Fit("mygaus_1_left","R");
mygaus_1_left->Draw("same");
auto l3 = new TLine(q[0]-1.5*iqr,0,q[0]-1.5*iqr,mygaus_1_left->Eval(q[0]-1.5*iqr));
l3->SetLineColor(kGreen); l3->SetLineWidth(2); l3->Draw("");
auto l1 = new TLine(q[0] ,0,q[0] ,mygaus_1_left->Eval(q[0]));
l1->SetLineWidth(2); l1->SetLineColor(kGreen); l1->Draw("");
h2->Fit("mygaus_1_right","R","");
mygaus_1_right->Draw("same");
auto l4 = new TLine(q[2]+1.5*iqr,0,q[2]+1.5*iqr,mygaus_1_left->Eval(q[2]+1.5*iqr));
l4->SetLineColor(kGreen); l4->SetLineWidth(2); l4->Draw("");
auto l5 = new TLine(q[2] ,0,q[2] ,mygaus_1_right->Eval(q[2]));
l5->SetLineWidth(2); l5->SetLineColor(kGreen); l5->Draw("");
h2->Fit("mygaus_1_middle","R");
mygaus_1_middle->Draw("same");
//In principal one could calculate these values by h2->Integral() as well
auto t = new TText(); t->SetTextFont(42);
t->DrawText(0,mygaus_1_middle->Eval(0)/2,"50%");
t->DrawText(-1.5,mygaus_1_middle->Eval(-1.5)/2,"24.65%");
t->DrawText(+1,mygaus_1_middle->Eval(+1.5)/2,"24.65%");
t->DrawText(q[0]-1.5*iqr,1000,Form("%.3f",q[0]-1.5*iqr))->SetTextAngle(90);
t->DrawText(q[2]+1.5*iqr,1000,Form("%.3f",q[2]+1.5*iqr))->SetTextAngle(90);
t->DrawText(q[0],1000,Form("%.3f",q[0]))->SetTextAngle(90);
t->DrawText(q[2],1000,Form("%.3f",q[2]))->SetTextAngle(90);
}
@ kGreen
Definition: Rtypes.h:64
float * q
Definition: THbookFile.cxx:87
char * Form(const char *fmt,...)
virtual void SetFillStyle(Style_t fstyle)
Set the fill area style.
Definition: TAttFill.h:39
virtual void SetLineWidth(Width_t lwidth)
Set the line width.
Definition: TAttLine.h:43
1-Dim function class
Definition: TF1.h:211
1-D histogram with an int per channel (see TH1 documentation)}
Definition: TH1.h:530
virtual Int_t GetQuantiles(Int_t nprobSum, Double_t *q, const Double_t *probSum=0)
Compute Quantiles for this histogram Quantile x_q of a probability distribution Function F is defined...
Definition: TH1.cxx:4451
TAxis * GetXaxis()
Get the behaviour adopted by the object about the statoverflows. See EStatOverflows for more informat...
Definition: TH1.h:316
virtual TFitResultPtr Fit(const char *formula, Option_t *option="", Option_t *goption="", Double_t xmin=0, Double_t xmax=0)
Fit histogram with function fname.
Definition: TH1.cxx:3808
TAxis * GetYaxis()
Definition: TH1.h:317
2-D histogram with an int per channel (see TH1 documentation)}
Definition: TH2.h:212
A simple line.
Definition: TLine.h:23
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition: TNamed.cxx:164
virtual void Draw(Option_t *option="")
Default Draw method for all objects.
Definition: TObject.cxx:195
This is the base class for the ROOT Random number generators.
Definition: TRandom.h:27
Base class for several text objects.
Definition: TText.h:23

The candle reduces the information coming from a whole distribution into few values. Independently from the number of entries or the significance of the underlying distribution a candle will always look like a candle. So candle plots should be used carefully in particular with unknown distributions. The definition of a candle is based on unbinned data. Here, candles are created from binned data. Because of this, the deviation is connected to the bin width used. The calculation of the quantiles normally done on unbinned data also. Because data are binned, this will only work the best possible way within the resolution of one bin

Because of all these facts one should take care that:

  • there are enough points per candle
  • the bin width is small enough (more bins will increase the maximum available resolution of the quantiles although there will be some bins with no entries)
  • never make a candle-plot if the underlying distribution is double-distributed
  • only create candles of distributions that are more-or-less gaussian (the MPV should be not too far away from the mean).

What a candle is made of

Since
ROOT version 6.07/05
The box

The box displays the position of the inter-quantile-range of the underlying distribution. The box contains 25% of the distribution below the median and 25% of the distribution above the median. If the underlying distribution is large enough and gaussian shaped the end-points of the box represent \( 0.6745\times\sigma \) (Where \( \sigma \) is the standard deviation of the gaussian). The width and the position of the box can be modified by SetBarWidth() and SetBarOffset(). The +-25% quantiles are calculated by the GetQuantiles() methods.

Since
ROOT version 6.11/01

Using the static function TCandle::SetBoxRange(double) the box definition will be overwritten. E.g. using a box range of 0.68 will redefine the area of the lower box edge to the upper box edge in order to cover 68% of the distribution illustrated by that candle. The static function will affect all candle-charts in the running program. Default is 0.5.

Using the static function TCandle::SetScaledCandle(bool) the width of the box (and the whole candle) can be influenced. Deactivated, the width is constant (to be set by SetBarWidth() ). Activated, the width of the boxes will be scaled to each other based on the amount of data in the corresponding candle, the maximum width can be influenced by SetBarWidth(). The static function will affect all candle-charts in the running program. Default is false. Scaling between multiple candle-charts (using "same" or THStack) is not supported, yet

The Median

For a sorted list of numbers, the median is the value in the middle of the list. E.g. if a sorted list is made of five numbers "1,2,3,6,7" 3 will be the median because it is in the middle of the list. If the number of entries is even the average of the two values in the middle will be used. As histograms are binned data, the situation is a bit more complex. The following example shows this:

void quantiles() {
auto h = new TH1I("h","h",10,0,10);
//h->Fill(3);
//h->Fill(3);
h->Fill(4);
h->Draw();
double p = 0.;
double q = 0.;
h->GetQuantiles(1,&q,&p);
cout << "Median is: " << q << std::endl;
}
#define h(i)
Definition: RSha256.hxx:106

Here the bin-width is 1.0. If the two Fill(3) are commented out, as there are currently, the example will return a calculated median of 4.5, because that's the bin center of the bin in which the value 4.0 has been dropped. If the two Fill(3) are not commented out, it will return 3.75, because the algorithm tries to evenly distribute the individual values of a bin with bin content > 0. It means the sorted list would be "3.25, 3.75, 4.5".

The consequence is a median of 3.75. This shows how important it is to use a small enough bin-width when using candle-plots on binned data. If the distribution is large enough and gaussian shaped the median will be exactly equal to the mean. The median can be shown as a line or as a circle or not shown at all.

In order to show the significance of the median notched candle plots apply a "notch" or narrowing of the box around the median. The significance is defined by \( 1.57\times\frac{iqr}{N} \) and will be represented as the size of the notch (where iqr is the size of the box and N is the number of entries of the whole distribution). Candle plots like these are usually called "notched candle plots".

In case the significance of the median is greater that the size of the box, the box will have an unnatural shape. Usually it means the chart has not enough data, or that representing this uncertainty is not useful

The Mean

The mean can be drawn as a dashed line or as a circle or not drawn at all. The mean is the arithmetic average of the values in the distribution. It is calculated using GetMean(). Because histograms are binned data, the mean value can differ from a calculation on the raw-data. If the distribution is large enough and gaussian shaped the mean will be exactly the median.

The Whiskers

The whiskers represent the part of the distribution not covered by the box. The upper 25% and the lower 25% of the distribution are located within the whiskers. Two representations are available.

  • A simple one (using w=1) defining the lower whisker from the lowest data value to the bottom of the box, and the upper whisker from the top of the box to the highest data value. In this representation the whisker-lines are dashed.
  • A more complex one having a further restriction. The whiskers are still connected to the box but their length cannot exceed \( 1.5\times iqr \). So it might be that the outermost part of the underlying distribution will not be covered by the whiskers. Usually these missing parts will be represented by the outliers (see points). Of course the upper and the lower whisker may differ in length. In this representation the whiskers are drawn as solid lines.
Since
ROOT version 6.11/01

Using the static function TCandle::SetWhiskerRange(double) the whisker definition w=1 will be overwritten. E.g. using a whisker-range of 0.95 and w=1 will redefine the area of the lower whisker to the upper whisker in order to cover 95% of the distribution inside that candle. The static function will affect all candle-charts in the running program. Default is 1.

If the distribution is large enough and gaussian shaped, the maximum length of the whisker will be located at \( \pm 2.698 \sigma \) (when using the 1.5*iqr-definition (w=2), where \( \sigma \) is the standard deviation (see picture above). In that case 99.3% of the total distribution will be covered by the box and the whiskers, whereas 0.7% are represented by the outliers.

The Anchors

The anchors have no special meaning in terms of statistical calculation. They mark the end of the whiskers and they have the width of the box. Both representation with and without anchors are common.

The Points

Depending on the configuration the points can have different meanings:

  • If p=1 the points represent the outliers. If they are shown, it means some parts of the underlying distribution are not covered by the whiskers. This can only occur when the whiskers are set to option w=2. Here the whiskers can have a maximum length of \( 1.5 \times iqr \). So any points outside the whiskers will be drawn as outliers. The outliers will be represented by crosses.
  • If p=2 all points in the distribution will be painted as crosses. This is useful for small datasets only (up to 10 or 20 points per candle). The outliers are shown along the candle. Because the underlying distribution is binned, is frequently occurs that a bin contains more than one value. Because of this the points will be randomly scattered within their bin along the candle axis. If the bin content for a bin is exactly 1 (usually this happens for the outliers) if will be drawn in the middle of the bin along the candle axis. As the maximum number of points per candle is limited by kNMax/2 on very large datasets scaling will be performed automatically. In that case one would loose all outliers because they have usually a bin content of 1 (and a bin content between 0 and 1 after the scaling). Because of this all bin contents between 0 and 1 - after the scaling - will be forced to be 1.
  • As the drawing of all values on large datasets can lead to big amounts of crosses, one can show all values as a scatter plot instead by choosing p=3. The points will be drawn as dots and will be scattered within the width of the candle. The color of the points will be the color of the candle-chart.
Other Options

Is is possible to combine all options of candle and violin plots with each other. E.g. a box-plot with a histogram.

How to use the candle-plots drawing option

There are six predefined candle-plot representations:

  • "CANDLEX1": Standard candle (whiskers cover the whole distribution)
  • "CANDLEX2": Standard candle with better whisker definition + outliers. It is a good compromise
  • "CANDLEX3": Like candle2 but with a mean as a circle. It is easier to distinguish mean and median
  • "CANDLEX4": Like candle3 but showing the uncertainty of the median as well (notched candle plots). For bigger datasets per candle
  • "CANDLEX5": Like candle2 but showing all data points. For very small datasets
  • "CANDLEX6": Like candle2 but showing all datapoints scattered. For huge datasets

The following picture shows how the six predefined representations look.

Example 1

Box and improved whisker, no mean, no median, no anchor no outliers

h1->Draw("CANDLEX(2001)");

Example 2

A Candle-definition like "CANDLEX2" (New standard candle with better whisker definition + outliers)

h1->Draw("CANDLEX(112111)");

Example 3

The following example shows how several candle plots can be super-imposed using the option SAME. Note that the bar-width and bar-offset are active on candle plots. Also the color, the line width, the size of the points and so on can be changed by the standard attribute setting methods such as SetLineColor() SetLineWidth().

void candleplot() {
TDatime dateBegin(2010,1,1,0,0,0);
TDatime dateEnd(2011,1,1,0,0,0);
auto h1 = new TH2I("h1","Machine A + B",12,dateBegin.Convert(),dateEnd.Convert(),1000,0,1000);
auto h2 = new TH2I("h2","Machine B",12,dateBegin.Convert(),dateEnd.Convert(),1000,0,1000);
h1->GetXaxis()->SetTimeFormat("%m/%y");
h1->GetXaxis()->SetTitle("Date [month/year]");
float Rand;
for (int i = dateBegin.Convert(); i < dateEnd.Convert(); i+=86400*30) {
for (int j = 0; j < 1000; j++) {
Rand = gRandom->Gaus(500+sin(i/10000000.)*100,50); h1->Fill(i,Rand);
Rand = gRandom->Gaus(500+sin(i/11000000.)*100,70); h2->Fill(i,Rand);
}
}
h1->SetBarWidth(0.4);
h1->SetBarOffset(-0.25);
h1->SetFillStyle(1001);
h2->SetBarWidth(0.4);
h2->SetBarOffset(0.25);
auto c1 = new TCanvas();
h1->Draw("candle2");
h2->Draw("candle3 same");
gPad->BuildLegend(0.78,0.695,0.980,0.935,"","f");
}
@ kYellow
Definition: Rtypes.h:64
double sin(double)
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition: TAttLine.h:40
virtual void SetTimeDisplay(Int_t value)
Definition: TAxis.h:161
virtual void SetTimeFormat(const char *format="")
Change the format used for time plotting.
Definition: TAxis.cxx:1002
This class stores the date and time with a precision of one second in an unsigned 32 bit word (950130...
Definition: TDatime.h:37
void SetTimeOffset(Double_t toffset)
Change the time offset for time plotting.
Definition: TStyle.cxx:1661

The VIOLIN option

A violin plot is a candle plot that also encodes the pdf information at each point.

Quartiles and mean are also represented at each point, with a marker and two lines.

In this implementation a TH2 is considered as a collection of TH1 along X (option VIOLIN or VIOLINX) or Y (option VIOLINY).

What a violin is made of

Since
ROOT version 6.09/02
The histogram

The histogram is typically drawn to both directions with respect to the middle-line of the corresponding bin. This can be achieved by using h=3. It is possible to draw a histogram only to one side (h=1, or h=2). The maximum number of bins in the histogram is limited to 500, if the number of bins in the used histogram is higher it will be rebinned automatically. The maximum height of the histogram can be modified by using SetBarWidth() and the position can be changed with SetBarOffset(). A solid fill style is recommended.

Since
ROOT version 6.11/01

Using the static function TCandle::SetScaledViolin(bool) the height of the histogram or the violin can be influenced. Activated, the height of the bins of the individual violins will be scaled with respect to each other, the maximum height can be influenced by SetBarWidth(). Deactivated, the height of the bin with the maximum content of each individual violin is set to a constant value using SetBarWidth(). The static function will affect all violin-charts in the running program. Default is true. Scaling between multiple violin-charts (using "same" or THStack) is not supported, yet.

The zero indicator line

Typical for violin charts is a line in the background over the whole histogram indicating the bins with zero entries. The zero indicator line can be activated with z=1. The line color will always be the same as the fill-color of the histogram.

The Mean

The Mean is illustrated with the same mechanism as used for candle plots. Usually a circle is used.

Whiskers

The whiskers are illustrated by the same mechanism as used for candle plots. There is only one difference. When using the simple whisker definition (w=1) and the zero indicator line (z=1), then the whiskers will be forced to be solid (usually hashed)

Points

The points are illustrated by the same mechanism as used for candle plots. E.g. VIOLIN2 uses better whisker definition (w=2) and outliers (p=1).

Other options

It is possible to combine all options of candle or violin plots with each other. E.g. a violin plot including a box-plot.

How to use the violin-plots drawing option

There are two predefined violin-plot representations:

  • "VIOLINX1": Standard violin (histogram, mean, whisker over full distribution, zero indicator line)
  • "VIOLINX2": Line VIOLINX1 both with better whisker definition + outliers.

A solid fill style is recommended for this plot (as opposed to a hollow or hashed style).

{
auto c1 = new TCanvas("c1","c1",600,400);
Int_t nx(6), ny(40);
double xmin(0.0), xmax(+6.0), ymin(0.0), ymax(+4.0);
auto hviolin = new TH2F("hviolin", "Option VIOLIN example", nx, xmin, xmax, ny, ymin, ymax);
TF1 f1("f1", "gaus", +0,0 +4.0);
double x,y;
for (Int_t iBin=1; iBin<hviolin->GetNbinsX(); ++iBin) {
double xc = hviolin->GetXaxis()->GetBinCenter(iBin);
f1.SetParameters(1, 2.0+TMath::Sin(1.0+xc), 0.2+0.1*(xc-xmin)/xmax);
for(Int_t i=0; i<10000; ++i){
x = xc;
y = f1.GetRandom();
hviolin->Fill(x, y);
}
}
hviolin->SetFillColor(kGray);
hviolin->SetMarkerStyle(20);
hviolin->SetMarkerSize(0.5);
hviolin->Draw("VIOLIN");
c1->Update();
}
@ kGray
Definition: Rtypes.h:63
float xmin
Definition: THbookFile.cxx:93
float ymin
Definition: THbookFile.cxx:93
float xmax
Definition: THbookFile.cxx:93
float ymax
Definition: THbookFile.cxx:93
virtual Double_t GetRandom()
Return a random number following this function shape.
Definition: TF1.cxx:2074
virtual void SetParameters(const Double_t *params)
Definition: TF1.h:638
TF1 * f1
Definition: legend1.C:11
Double_t Sin(Double_t)
Definition: TMath.h:627

The next example illustrates a time development of a certain value:

void candledecay()
{
auto c1 = new TCanvas("c1","Candle Decay",800,600);
c1->Divide(2,1);
auto rng = new TRandom();
auto h1 = new TH2I("h1","Decay",1000,0,1000,20,0,20);
float myRand;
for (int i = 0; i < 19; i++) {
for (int j = 0; j < 1000000; j++) {
myRand = rng->Gaus(350+i*8,20+2*i);
h1->Fill(myRand,i);
}
}
h1->GetYaxis()->SetTitle("time");
h1->GetXaxis()->SetTitle("probability density");
c1->cd(1);
h1->Draw("violiny(112000000)");
c1->cd(2);
auto h2 = (TH2I*)h1->Clone("h2");
h2->SetBarWidth(0.8);
h2->DrawCopy("candley2");
}
@ kBlue
Definition: Rtypes.h:64

The TEXT and TEXTnn Option

For each bin the content is printed. The text attributes are:

  • text font = current TStyle font (gStyle->SetTextFont()).
  • text size = 0.02*padheight*markersize (if h is the histogram drawn with the option TEXT the marker size can be changed with h->SetMarkerSize(markersize)).
  • text color = marker color.

By default the format g is used. This format can be redefined by calling gStyle->SetPaintTextFormat().

It is also possible to use TEXTnn in order to draw the text with the angle nn (0 < nn < 90).

For 2D histograms the text is plotted in the center of each non empty cells. It is possible to plot empty cells by calling gStyle->SetHistMinimumZero() or providing MIN0 draw option. For 1D histogram the text is plotted at a y position equal to the bin content.

For 2D histograms when the option "E" (errors) is combined with the option text ("TEXTE"), the error for each bin is also printed.

{
auto c01 = new TCanvas("c01","c01",700,400);
c01->Divide(2,1);
auto htext1 = new TH1F("htext1","Option TEXT on 1D histograms ",10,-4,4);
auto htext2 = new TH2F("htext2","Option TEXT on 2D histograms ",10,-4,4,10,-20,20);
float px, py;
for (Int_t i = 0; i < 25000; i++) {
gRandom->Rannor(px,py);
htext1->Fill(px,0.1);
htext2->Fill(px,5*py,0.1);
}
htext2->SetMarkerSize(1.8);
c01->cd(1);
htext2->Draw("TEXT45");
c01->cd(2);
htext1->Draw();
htext1->Draw("HIST TEXT0 SAME");
}
void SetPaintTextFormat(const char *format="g")
Definition: TStyle.h:364
Since
ROOT version 6.07/07:

In case several histograms are drawn on top ot each other (using option SAME), the text can be shifted using SetBarOffset(). It specifies an offset for the text position in each cell, in percentage of the bin width.

{
auto c03 = new TCanvas("c03","c03",700,400);
auto htext3 = new TH2F("htext3","Several 2D histograms drawn with option TEXT",10,-4,4,10,-20,20);
auto htext4 = new TH2F("htext4","htext4",10,-4,4,10,-20,20);
auto htext5 = new TH2F("htext5","htext5",10,-4,4,10,-20,20);
float px, py;
for (Int_t i = 0; i < 25000; i++) {
gRandom->Rannor(px,py);
htext3->Fill(4*px,20*py,0.1);
htext4->Fill(4*px,20*py,0.5);
htext5->Fill(4*px,20*py,1.0);
}
htext4->SetMarkerSize(1.8);
htext5->SetMarkerSize(1.8);
htext5->SetMarkerColor(kRed);
htext3->Draw("COL");
htext4->SetBarOffset(0.2);
htext4->Draw("TEXT SAME");
htext5->SetBarOffset(-0.2);
htext5->Draw("TEXT SAME");
}

In the case of profile histograms it is possible to print the number of entries instead of the bin content. It is enough to combine the option "E" (for entries) with the option "TEXT".

{
auto c02 = new TCanvas("c02","c02",700,400);
c02->Divide(2,1);
auto profile = new TProfile("profile","profile",10,0,10);
profile->SetMarkerSize(2.2);
profile->Fill(0.5,1);
profile->Fill(1.5,2);
profile->Fill(2.5,3);
profile->Fill(3.5,4);
profile->Fill(4.5,5);
profile->Fill(5.5,5);
profile->Fill(6.5,4);
profile->Fill(7.5,3);
profile->Fill(8.5,2);
profile->Fill(9.5,1);
c02->cd(1); profile->Draw("HIST TEXT0");
c02->cd(2); profile->Draw("HIST TEXT0E");
}
Profile Histogram.
Definition: TProfile.h:32

The CONTour options

The following contour options are supported:

Option Description
"CONT" Draw a contour plot (same as CONT0).
"CONT0" Draw a contour plot using surface colors to distinguish contours.
"CONT1" Draw a contour plot using the line colors to distinguish contours.
"CONT2" Draw a contour plot using the line styles to distinguish contours.
"CONT3" Draw a contour plot solid lines for all contours.
"CONT4" Draw a contour plot using surface colors (SURF option at theta = 0).
"CONT5" Draw a contour plot using Delaunay triangles.

The following example shows a 2D histogram plotted with the option CONTZ. The option CONT draws a contour plot using surface colors to distinguish contours. Combined with the option CONT (or CONT0), the option Z allows to display the color palette defined by gStyle->SetPalette().

{
auto c1 = new TCanvas("c1","c1",600,400);
auto hcontz = new TH2F("hcontz","Option CONTZ example ",40,-4,4,40,-20,20);
float px, py;
for (Int_t i = 0; i < 25000; i++) {
gRandom->Rannor(px,py);
hcontz->Fill(px-1,5*py);
hcontz->Fill(2+0.5*px,2*py-10.,0.1);
}
hcontz->Draw("CONTZ");
}

The following example shows a 2D histogram plotted with the option CONT1Z. The option CONT1 draws a contour plot using the line colors to distinguish contours. Combined with the option CONT1, the option Z allows to display the color palette defined by gStyle->SetPalette().

{
auto c1 = new TCanvas("c1","c1",600,400);
auto hcont1 = new TH2F("hcont1","Option CONT1Z example ",40,-4,4,40,-20,20);
float px, py;
for (Int_t i = 0; i < 25000; i++) {
gRandom->Rannor(px,py);
hcont1->Fill(px-1,5*py);
hcont1->Fill(2+0.5*px,2*py-10.,0.1);
}
hcont1->Draw("CONT1Z");
}

The following example shows a 2D histogram plotted with the option CONT2. The option CONT2 draws a contour plot using the line styles to distinguish contours.

{
auto c1 = new TCanvas("c1","c1",600,400);
auto hcont2 = new TH2F("hcont2","Option CONT2 example ",40,-4,4,40,-20,20);
float px, py;
for (Int_t i = 0; i < 25000; i++) {
gRandom->Rannor(px,py);
hcont2->Fill(px-1,5*py);
hcont2->Fill(2+0.5*px,2*py-10.,0.1);
}
hcont2->Draw("CONT2");
}

The following example shows a 2D histogram plotted with the option CONT3. The option CONT3 draws contour plot solid lines for all contours.

{
auto c1 = new TCanvas("c1","c1",600,400);
auto hcont3 = new TH2F("hcont3","Option CONT3 example ",40,-4,4,40,-20,20);
float px, py;
for (Int_t i = 0; i < 25000; i++) {
gRandom->Rannor(px,py);
hcont3->Fill(px-1,5*py);
hcont3->Fill(2+0.5*px,2*py-10.,0.1);
}
hcont3->Draw("CONT3");
}

The following example shows a 2D histogram plotted with the option CONT4. The option CONT4 draws a contour plot using surface colors to distinguish contours (SURF option at theta = 0). Combined with the option CONT (or CONT0), the option Z allows to display the color palette defined by gStyle->SetPalette().

{
auto c1 = new TCanvas("c1","c1",600,400);
auto hcont4 = new TH2F("hcont4","Option CONT4Z example ",40,-4,4,40,-20,20);
float px, py;
for (Int_t i = 0; i < 25000; i++) {
gRandom->Rannor(px,py);
hcont4->Fill(px-1,5*py);
hcont4->Fill(2+0.5*px,2*py-10.,0.1);
}
hcont4->Draw("CONT4Z");
}

The default number of contour levels is 20 equidistant levels and can be changed with TH1::SetContour() or TStyle::SetNumberContours().

The LIST option

When option LIST is specified together with option CONT, the points used to draw the contours are saved in TGraph objects:

h->Draw("CONT LIST");
gPad->Update();

The contour are saved in TGraph objects once the pad is painted. Therefore to use this functionality in a macro, gPad->Update() should be performed after the histogram drawing. Once the list is built, the contours are accessible in the following way:

TObjArray *contours = (TObjArray*)gROOT->GetListOfSpecials()->FindObject("contours");
Int_t ncontours     = contours->GetSize();
TList *list         = (TList*)contours->At(i);

Where i is a contour number, and list contains a list of TGraph objects. For one given contour, more than one disjoint polyline may be generated. The number of TGraphs per contour is given by:

list->GetSize();

To access the first graph in the list one should do:

TGraph *gr1 = (TGraph*)list->First();

The following example (ContourList.C) shows how to use this functionality.

Double_t SawTooth(Double_t x, Double_t WaveLen);
TCanvas *ContourList(){
const Double_t PI = TMath::Pi();
TCanvas* c = new TCanvas("c","Contour List",0,0,600,600);
c->SetRightMargin(0.15);
c->SetTopMargin(0.15);
Int_t i, j;
Int_t nZsamples = 80;
Int_t nPhiSamples = 80;
Double_t HofZwavelength = 4.0; // 4 meters
Double_t dZ = HofZwavelength/(Double_t)(nZsamples - 1);
Double_t dPhi = 2*PI/(Double_t)(nPhiSamples - 1);
TArrayD z(nZsamples);
TArrayD HofZ(nZsamples);
TArrayD phi(nPhiSamples);
TArrayD FofPhi(nPhiSamples);
// Discretized Z and Phi Values
for ( i = 0; i < nZsamples; i++) {
z[i] = (i)*dZ - HofZwavelength/2.0;
HofZ[i] = SawTooth(z[i], HofZwavelength);
}
for(Int_t i=0; i < nPhiSamples; i++){
phi[i] = (i)*dPhi;
FofPhi[i] = sin(phi[i]);
}
// Create Histogram
TH2D *HistStreamFn = new TH2D("HstreamFn",
"#splitline{Histogram with negative and positive contents. Six contours are defined.}{It is plotted with options CONT LIST to retrieve the contours points in TGraphs}",
nZsamples, z[0], z[nZsamples-1], nPhiSamples, phi[0], phi[nPhiSamples-1]);
// Load Histogram Data
for (Int_t i = 0; i < nZsamples; i++) {
for(Int_t j = 0; j < nPhiSamples; j++){
HistStreamFn->SetBinContent(i,j, HofZ[i]*FofPhi[j]);
}
}
gStyle->SetTitleW(0.99);
gStyle->SetTitleH(0.08);
Double_t contours[6];
contours[0] = -0.7;
contours[1] = -0.5;
contours[2] = -0.1;
contours[3] = 0.1;
contours[4] = 0.4;
contours[5] = 0.8;
HistStreamFn->SetContour(6, contours);
// Draw contours as filled regions, and Save points
HistStreamFn->Draw("CONT Z LIST");
c->Update(); // Needed to force the plotting and retrieve the contours in TGraphs
// Get Contours
TObjArray *conts = (TObjArray*)gROOT->GetListOfSpecials()->FindObject("contours");
TList* contLevel = NULL;
TGraph* curv = NULL;
TGraph* gc = NULL;
Int_t nGraphs = 0;
Int_t TotalConts = 0;
if (conts == NULL){
printf("*** No Contours Were Extracted!\n");
TotalConts = 0;
return 0;
} else {
TotalConts = conts->GetSize();
}
printf("TotalConts = %d\n", TotalConts);
for(i = 0; i < TotalConts; i++){
contLevel = (TList*)conts->At(i);
printf("Contour %d has %d Graphs\n", i, contLevel->GetSize());
nGraphs += contLevel->GetSize();
}
nGraphs = 0;
TCanvas* c1 = new TCanvas("c1","Contour List",610,0,600,600);
c1->SetTopMargin(0.15);
TH2F *hr = new TH2F("hr",
"#splitline{Negative contours are returned first (highest to lowest). Positive contours are returned from}{lowest to highest. On this plot Negative contours are drawn in red and positive contours in blue.}",
2, -2, 2, 2, 0, 6.5);
hr->Draw();
Double_t xval0, yval0, zval0;
l.SetTextSize(0.03);
char val[20];
for(i = 0; i < TotalConts; i++){
contLevel = (TList*)conts->At(i);
if (i<3) zval0 = contours[2-i];
else zval0 = contours[i];
printf("Z-Level Passed in as: Z = %f\n", zval0);
// Get first graph from list on curves on this level
curv = (TGraph*)contLevel->First();
for(j = 0; j < contLevel->GetSize(); j++){
curv->GetPoint(0, xval0, yval0);
if (zval0<0) curv->SetLineColor(kRed);
if (zval0>0) curv->SetLineColor(kBlue);
nGraphs ++;
printf("\tGraph: %d -- %d Elements\n", nGraphs,curv->GetN());
// Draw clones of the graphs to avoid deletions in case the 1st
// pad is redrawn.
gc = (TGraph*)curv->Clone();
gc->Draw("C");
sprintf(val,"%g",zval0);
l.DrawLatex(xval0,yval0,val);
curv = (TGraph*)contLevel->After(curv); // Get Next graph
}
}
c1->Update();
printf("\n\n\tExtracted %d Contours and %d Graphs \n", TotalConts, nGraphs );
return c1;
}
Double_t SawTooth(Double_t x, Double_t WaveLen){
// This function is specific to a sawtooth function with period
// WaveLen, symmetric about x = 0, and with amplitude = 1. Each segment
// is 1/4 of the wavelength.
//
// |
// /\ |
// / \ |
// / \ |
// / \
// /--------\--------/------------
// |\ /
// | \ /
// | \ /
// | \/
//
if ( (x < -WaveLen/2) || (x > WaveLen/2)) y = -99999999; // Error X out of bounds
if (x <= -WaveLen/4) {
y = x + 2.0;
} else if ((x > -WaveLen/4) && (x <= WaveLen/4)) {
y = -x ;
} else if (( x > WaveLen/4) && (x <= WaveLen/2)) {
y = x - 2.0;
}
return y;
}
#define PI
Array of doubles (64 bits per element).
Definition: TArrayD.h:27
virtual Int_t GetSize() const
Return the capacity of the collection, i.e.
Definition: TCollection.h:182
A Graph is a graphics object made of two arrays X and Y with npoints each.
Definition: TGraph.h:41
Int_t GetN() const
Definition: TGraph.h:123
virtual void Draw(Option_t *chopt="")
Draw this graph with its current attributes.
Definition: TGraph.cxx:753
virtual Int_t GetPoint(Int_t i, Double_t &x, Double_t &y) const
Get x and y values for point number i.
Definition: TGraph.cxx:1586
virtual void SetContour(Int_t nlevels, const Double_t *levels=0)
Set the number and values of contour levels.
Definition: TH1.cxx:7935
2-D histogram with a double per channel (see TH1 documentation)}
Definition: TH2.h:292
virtual void SetBinContent(Int_t bin, Double_t content)
Set bin content.
Definition: TH2.cxx:2441
A doubly linked list.
Definition: TList.h:44
virtual TObject * After(const TObject *obj) const
Returns the object after object obj.
Definition: TList.cxx:327
virtual TObject * First() const
Return the first object in the list. Returns 0 when list is empty.
Definition: TList.cxx:656
virtual TObject * Clone(const char *newname="") const
Make a clone of an object using the Streamer facility.
Definition: TNamed.cxx:74
An array of TObjects.
Definition: TObjArray.h:37
TObject * At(Int_t idx) const
Definition: TObjArray.h:166
void SetTitleW(Float_t w=0)
Definition: TStyle.h:393
void SetTitleH(Float_t h=0)
Definition: TStyle.h:394
constexpr Double_t Pi()
Definition: TMath.h:38
auto * l
Definition: textangle.C:4

The AITOFF, MERCATOR, SINUSOIDAL and PARABOLIC options

The following options select the CONT4 option and are useful for sky maps or exposure maps (earth.C).

Option Description
"AITOFF" Draw a contour via an AITOFF projection.
"MERCATOR" Draw a contour via an Mercator projection.
"SINUSOIDAL" Draw a contour via an Sinusoidal projection.
"PARABOLIC" Draw a contour via an Parabolic projection.
TCanvas *earth(){
TCanvas *c1 = new TCanvas("c1","earth_projections",700,700);
c1->Divide(2,2);
TH2F *ha = new TH2F("ha","Aitoff", 180, -180, 180, 179, -89.5, 89.5);
TH2F *hm = new TH2F("hm","Mercator", 180, -180, 180, 161, -80.5, 80.5);
TH2F *hs = new TH2F("hs","Sinusoidal",180, -180, 180, 181, -90.5, 90.5);
TH2F *hp = new TH2F("hp","Parabolic", 180, -180, 180, 181, -90.5, 90.5);
TString dat = gROOT->GetTutorialDir();
dat.Append("/graphics/earth.dat");
dat.ReplaceAll("/./","/");
ifstream in;
in.open(dat.Data());
while (1) {
in >> x >> y;
if (!in.good()) break;
ha->Fill(x,y, 1);
hm->Fill(x,y, 1);
hs->Fill(x,y, 1);
hp->Fill(x,y, 1);
}
in.close();
c1->cd(1); ha->Draw("aitoff");
c1->cd(2); hm->Draw("mercator");
c1->cd(3); hs->Draw("sinusoidal");
c1->cd(4); hp->Draw("parabolic");
return c1;
}
float Float_t
Definition: RtypesCore.h:53
Int_t Fill(Double_t)
Invalid Fill method.
Definition: TH2.cxx:292
const char * Data() const
Definition: TString.h:364
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition: TString.h:687
TString & Append(const char *cs)
Definition: TString.h:559

The LEGO options

In a lego plot the cell contents are drawn as 3-d boxes. The height of each box is proportional to the cell content. The lego aspect is control with the following options:

Option Description
"LEGO" Draw a lego plot using the hidden lines removal technique.
"LEGO1" Draw a lego plot using the hidden surface removal technique.
"LEGO2" Draw a lego plot using colors to show the cell contents.
"LEGO3" Draw a lego plot with hidden surface removal, like LEGO1 but the border lines of each lego-bar are not drawn.
"LEGO4" Draw a lego plot with hidden surface removal, like LEGO1 but without the shadow effect on each lego-bar.
"0" When used with any LEGO option, the empty bins are not drawn.

See the limitations with the option "SAME".

Line attributes can be used in lego plots to change the edges' style.

The following example shows a 2D histogram plotted with the option LEGO. The option LEGO draws a lego plot using the hidden lines removal technique.

{
auto c2 = new TCanvas("c2","c2",600,400);
auto hlego = new TH2F("hlego","Option LEGO example ",40,-4,4,40,-20,20);
float px, py;
for (Int_t i = 0; i < 25000; i++) {
gRandom->Rannor(px,py);
hlego->Fill(px-1,5*py);
hlego->Fill(2+0.5*px,2*py-10.,0.1);
}
hlego->Draw("LEGO");
}
return c2
Definition: legend2.C:14

The following example shows a 2D histogram plotted with the option LEGO1. The option LEGO1 draws a lego plot using the hidden surface removal technique. Combined with any LEGOn option, the option 0 allows to not drawn the empty bins.

{
auto c2 = new TCanvas("c2","c2",600,400);
auto hlego1 = new TH2F("hlego1","Option LEGO1 example (with option 0) ",40,-4,4,40,-20,20);
float px, py;
for (Int_t i = 0; i < 25000; i++) {
gRandom->Rannor(px,py);
hlego1->Fill(px-1,5*py);
hlego1->Fill(2+0.5*px,2*py-10.,0.1);
}
hlego1->SetFillColor(kYellow);
hlego1->Draw("LEGO1 0");
}

The following example shows a 2D histogram plotted with the option LEGO3. Like the option LEGO1, the option LEGO3 draws a lego plot using the hidden surface removal technique but doesn't draw the border lines of each individual lego-bar. This is very useful for histograms having many bins. With such histograms the option LEGO1 gives a black image because of the border lines. This option also works with stacked legos.

{
auto c2 = new TCanvas("c2","c2",600,400);
auto hlego3 = new TH2F("hlego3","Option LEGO3 example",40,-4,4,40,-20,20);
float px, py;
for (Int_t i = 0; i < 25000; i++) {
gRandom->Rannor(px,py);
hlego3->Fill(px-1,5*py);
hlego3->Fill(2+0.5*px,2*py-10.,0.1);
}
hlego3->SetFillColor(kRed);
hlego3->Draw("LEGO3");
}

The following example shows a 2D histogram plotted with the option LEGO2. The option LEGO2 draws a lego plot using colors to show the cell contents. Combined with the option LEGO2, the option Z allows to display the color palette defined by gStyle->SetPalette().

{
auto c2 = new TCanvas("c2","c2",600,400);
auto hlego2 = new TH2F("hlego2","Option LEGO2Z example ",40,-4,4,40,-20,20);
float px, py;
for (Int_t i = 0; i < 25000; i++) {
gRandom->Rannor(px,py);
hlego2->Fill(px-1,5*py);
hlego2->Fill(2+0.5*px,2*py-10.,0.1);
}
hlego2->Draw("LEGO2Z");
}

The "SURFace" options

In a surface plot, cell contents are represented as a mesh. The height of the mesh is proportional to the cell content.

Option Description
"SURF" Draw a surface plot using the hidden line removal technique.
"SURF1" Draw a surface plot using the hidden surface removal technique.
"SURF2" Draw a surface plot using colors to show the cell contents.
"SURF3" Same as SURF with an additional filled contour plot on top.
"SURF4" Draw a surface using the Gouraud shading technique.
"SURF5" Used with one of the options CYL, PSR and CYL this option allows to draw a a filled contour plot.
"SURF6" This option should not be used directly. It is used internally when the CONT is used with option the option SAME on a 3D plot.
"SURF7" Same as SURF2 with an additional line contour plot on top.

See the limitations with the option "SAME".

The following example shows a 2D histogram plotted with the option SURF. The option SURF draws a lego plot using the hidden lines removal technique.

{
auto c2 = new TCanvas("c2","c2",600,400);
auto hsurf = new TH2F("hsurf","Option SURF example ",30,-4,4,30,-20,20);
float px, py;
for (Int_t i = 0; i < 25000; i++) {
gRandom->Rannor(px,py);
hsurf->Fill(px-1,5*py);
hsurf->Fill(2+0.5*px,2*py-10.,0.1);
}
hsurf->Draw("SURF");
}

The following example shows a 2D histogram plotted with the option SURF1. The option SURF1 draws a surface plot using the hidden surface removal technique. Combined with the option SURF1, the option Z allows to display the color palette defined by gStyle->SetPalette().

{
auto c2 = new TCanvas("c2","c2",600,400);
auto hsurf1 = new TH2F("hsurf1","Option SURF1 example ",30,-4,4,30,-20,20);
float px, py;
for (Int_t i = 0; i < 25000; i++) {
gRandom->Rannor(px,py);
hsurf1->Fill(px-1,5*py);
hsurf1->Fill(2+0.5*px,2*py-10.,0.1);
}
hsurf1->Draw("SURF1");
}

The following example shows a 2D histogram plotted with the option SURF2. The option SURF2 draws a surface plot using colors to show the cell contents. Combined with the option SURF2, the option Z allows to display the color palette defined by gStyle->SetPalette().

{
auto c2 = new TCanvas("c2","c2",600,400);
auto hsurf2 = new TH2F("hsurf2","Option SURF2 example ",30,-4,4,30,-20,20);
float px, py;
for (Int_t i = 0; i < 25000; i++) {
gRandom->Rannor(px,py);
hsurf2->Fill(px-1,5*py);
hsurf2->Fill(2+0.5*px,2*py-10.,0.1);
}
hsurf2->Draw("SURF2");
}

The following example shows a 2D histogram plotted with the option SURF3. The option SURF3 draws a surface plot using the hidden line removal technique with, in addition, a filled contour view drawn on the top. Combined with the option SURF3, the option Z allows to display the color palette defined by gStyle->SetPalette().

{
auto c2 = new TCanvas("c2","c2",600,400);
auto hsurf3 = new TH2F("hsurf3","Option SURF3 example ",30,-4,4,30,-20,20);
float px, py;
for (Int_t i = 0; i < 25000; i++) {
gRandom->Rannor(px,py);
hsurf3->Fill(px-1,5*py);
hsurf3->Fill(2+0.5*px,2*py-10.,0.1);
}
hsurf3->Draw("SURF3");
}

The following example shows a 2D histogram plotted with the option SURF4. The option SURF4 draws a surface using the Gouraud shading technique.

{
auto c2 = new TCanvas("c2","c2",600,400);
auto hsurf4 = new TH2F("hsurf4","Option SURF4 example ",30,-4,4,30,-20,20);
float px, py;
for (Int_t i = 0; i < 25000; i++) {
gRandom->Rannor(px,py);
hsurf4->Fill(px-1,5*py);
hsurf4->Fill(2+0.5*px,2*py-10.,0.1);
}
hsurf4->SetFillColor(kOrange);
hsurf4->Draw("SURF4");
}
@ kOrange
Definition: Rtypes.h:65

The following example shows a 2D histogram plotted with the option SURF5 CYL. Combined with the option SURF5, the option Z allows to display the color palette defined by gStyle->SetPalette().

{
auto c2 = new TCanvas("c2","c2",600,400);
auto hsurf5 = new TH2F("hsurf4","Option SURF5 example ",30,-4,4,30,-20,20);
float px, py;
for (Int_t i = 0; i < 25000; i++) {
gRandom->Rannor(px,py);
hsurf5->Fill(px-1,5*py);
hsurf5->Fill(2+0.5*px,2*py-10.,0.1);
}
hsurf5->Draw("SURF5 CYL");
}

The following example shows a 2D histogram plotted with the option SURF7. The option SURF7 draws a surface plot using the hidden surfaces removal technique with, in addition, a line contour view drawn on the top. Combined with the option SURF7, the option Z allows to display the color palette defined by gStyle->SetPalette().

{
auto c2 = new TCanvas("c2","c2",600,400);
auto hsurf7 = new TH2F("hsurf3","Option SURF7 example ",30,-4,4,30,-20,20);
float px, py;
for (Int_t i = 0; i < 25000; i++) {
gRandom->Rannor(px,py);
hsurf7->Fill(px-1,5*py);
hsurf7->Fill(2+0.5*px,2*py-10.,0.1);
}
hsurf7->Draw("SURF7");
}

As shown in the following example, when a contour plot is painted on top of a surface plot using the option SAME, the contours appear in 3D on the surface.

{
auto c20=new TCanvas("c20","c20",600,400);
int NBins = 50;
double d = 2;
auto hsc = new TH2F("hsc", "Surface and contour with option SAME ", NBins, -d, d, NBins, -d, d);
for (int bx = 1; bx <= NBins; ++bx) {
for (int by = 1; by <= NBins; ++by) {
double x = hsc->GetXaxis()->GetBinCenter(bx);
double y = hsc->GetYaxis()->GetBinCenter(by);
hsc->SetBinContent(bx, by, exp(-x*x)*exp(-y*y));
}
}
hsc->Draw("surf2");
hsc->Draw("CONT1 SAME");
}
#define d(i)
Definition: RSha256.hxx:102
double exp(double)

Cylindrical, Polar, Spherical and PseudoRapidity/Phi options

Legos and surfaces plots are represented by default in Cartesian coordinates. Combined with any LEGOn or SURFn options the following options allow to draw a lego or a surface in other coordinates systems.

Option Description
"CYL" Use Cylindrical coordinates. The X coordinate is mapped on the angle and the Y coordinate on the cylinder length.
"POL" Use Polar coordinates. The X coordinate is mapped on the angle and the Y coordinate on the radius.
"SPH" Use Spherical coordinates. The X coordinate is mapped on the latitude and the Y coordinate on the longitude.
"PSR" Use PseudoRapidity/Phi coordinates. The X coordinate is mapped on Phi.

WARNING: Axis are not drawn with these options.

The following example shows the same histogram as a lego plot is the four different coordinates systems.

{
auto c3 = new TCanvas("c3","c3",600,400);
c3->Divide(2,2);
auto hlcc = new TH2F("hlcc","Cylindrical coordinates",20,-4,4,20,-20,20);
float px, py;
for (Int_t i = 0; i < 25000; i++) {
gRandom->Rannor(px,py);
hlcc->Fill(px-1,5*py);
hlcc->Fill(2+0.5*px,2*py-10.,0.1);
}
hlcc->SetFillColor(kYellow);
c3->cd(1); hlcc->Draw("LEGO1 CYL");
c3->cd(2); auto hlpc = (TH2F*) hlcc->DrawClone("LEGO1 POL");
hlpc->SetTitle("Polar coordinates");
c3->cd(3); auto hlsc = (TH2F*) hlcc->DrawClone("LEGO1 SPH");
hlsc->SetTitle("Spherical coordinates");
c3->cd(4); auto hlprpc = (TH2F*) hlcc->DrawClone("LEGO1 PSR");
hlprpc->SetTitle("PseudoRapidity/Phi coordinates");
}
return c3
Definition: legend3.C:15

The following example shows the same histogram as a surface plot is the four different coordinates systems.

{
auto c4 = new TCanvas("c4","c4",600,400);
c4->Divide(2,2);
auto hscc = new TH2F("hscc","Cylindrical coordinates",20,-4,4,20,-20,20);
float px, py;
for (Int_t i = 0; i < 25000; i++) {
gRandom->Rannor(px,py);
hscc->Fill(px-1,5*py);
hscc->Fill(2+0.5*px,2*py-10.,0.1);
}
c4->cd(1); hscc->Draw("SURF1 CYL");
c4->cd(2); auto hspc = (TH2F*) hscc->DrawClone("SURF1 POL");
hspc->SetTitle("Polar coordinates");
c4->cd(3); auto hssc = (TH2F*) hscc->DrawClone("SURF1 SPH");
hssc->SetTitle("Spherical coordinates");
c4->cd(4); auto hsprpc = (TH2F*) hscc->DrawClone("SURF1 PSR");
hsprpc->SetTitle("PseudoRapidity/Phi coordinates");
}

Base line for bar-charts and lego plots

By default the base line used to draw the boxes for bar-charts and lego plots is the histogram minimum. It is possible to force this base line to be 0, using MIN0 draw option or with the command:

gStyle->SetHistMinimumZero();
{
auto c5 = new TCanvas("c5","c5",700,400);
c5->Divide(2,1);
auto hz1 = new TH1F("hz1","Bar-chart drawn from 0",20,-3,3);
auto hz2 = new TH2F("hz2","Lego plot drawn from 0",20,-3,3,20,-3,3);
Int_t i;
double x,y;
hz1->SetFillColor(kBlue);
hz2->SetFillColor(kBlue);
for (i=0;i<10000;i++) {
x = gRandom->Gaus(0,1);
y = gRandom->Gaus(0,1);
if (x>0) {
hz1->Fill(x,1);
hz2->Fill(x,y,1);
} else {
hz1->Fill(x,-1);
hz2->Fill(x,y,-2);
}
}
c5->cd(1); hz1->Draw("bar2 min0");
c5->cd(2); hz2->Draw("lego1 min0");
}

This option also works for horizontal plots. The example given in the section "The bar chart option" appears as follow:

{
int i;
const Int_t nx = 8;
string os_X[nx] = {"8","32","128","512","2048","8192","32768","131072"};
float d_35_0[nx] = {0.75, -3.30, -0.92, 0.10, 0.08, -1.69, -1.29, -2.37};
float d_35_1[nx] = {1.01, -3.02, -0.65, 0.37, 0.34, -1.42, -1.02, -2.10};
auto cbh = new TCanvas("cbh","cbh",400,600);
cbh->SetGrid();
auto h1bh = new TH1F("h1bh","Option HBAR centered on 0",nx,0,nx);
h1bh->SetFillColor(4);
h1bh->SetBarWidth(0.4);
h1bh->SetBarOffset(0.1);
h1bh->SetStats(0);
h1bh->SetMinimum(-5);
h1bh->SetMaximum(5);
for (i=1; i<=nx; i++) {
h1bh->Fill(os_X[i-1].c_str(), d_35_0[i-1]);
h1bh->GetXaxis()->SetBinLabel(i,os_X[i-1].c_str());
}
h1bh->Draw("hbar min0");
auto h2bh = new TH1F("h2bh","h2bh",nx,0,nx);
h2bh->SetFillColor(38);
h2bh->SetBarWidth(0.4);
h2bh->SetBarOffset(0.5);
h2bh->SetStats(0);
for (i=1;i<=nx;i++) h2bh->Fill(os_X[i-1].c_str(), d_35_1[i-1]);
h2bh->Draw("hbar min0 same");
}

TH2Poly Drawing

The following options are supported:

Option Description
"SCAT" Draw a scatter plot (default).
"COL" Draw a color plot. All the bins are painted even the empty bins.
"COLZ" Same as "COL". In addition the color palette is also drawn.
"0" When used with any COL options, the empty bins are not drawn.
"TEXT" Draw bin contents as text (format set via gStyle->SetPaintTextFormat).
"TEXTN" Draw bin names as text.
"TEXTnn" Draw bin contents as text at angle nn (0 < nn < 90).
"L" Draw the bins boundaries as lines. The lines attributes are the TGraphs ones.
"P" Draw the bins boundaries as markers. The markers attributes are the TGraphs ones.
"F" Draw the bins boundaries as filled polygons. The filled polygons attributes are the TGraphs ones.

TH2Poly can be drawn as a color plot (option COL). TH2Poly bins can have any shapes. The bins are defined as graphs. The following macro is a very simple example showing how to book a TH2Poly and draw it.

{
auto ch2p1 = new TCanvas("ch2p1","ch2p1",600,400);
auto h2p = new TH2Poly();
h2p->SetName("h2poly_name");
h2p->SetTitle("h2poly_title");
double px1[] = {0, 5, 6};
double py1[] = {0, 0, 5};
double px2[] = {0, -1, -1, 0};
double py2[] = {0, 0, -1, 3};
double px3[] = {4, 3, 0, 1, 2.4};
double py3[] = {4, 3.7, 1, 3.7, 2.5};
h2p->AddBin(3, px1, py1);
h2p->AddBin(4, px2, py2);
h2p->AddBin(5, px3, py3);
h2p->Fill(0.1, 0.01, 3);
h2p->Fill(-0.5, -0.5, 7);
h2p->Fill(-0.7, -0.5, 1);
h2p->Fill(1, 3, 1.5);
double fx[] = {0.1, -0.5, -0.7, 1};
double fy[] = {0.01, -0.5, -0.5, 3};
double fw[] = {3, 1, 1, 1.5};
h2p->FillN(4, fx, fy, fw);
h2p->Draw("col");
}
2D Histogram with Polygonal Bins
Definition: TH2Poly.h:66

Rectangular bins are a frequent case. The special version of the AddBin method allows to define them more easily like shown in the following example (th2polyBoxes.C).

TCanvas *th2polyBoxes() {
TCanvas *ch2p2 = new TCanvas("ch2p2","ch2p2",600,400);
TH2Poly *h2p = new TH2Poly();
h2p->SetName("Boxes");
h2p->SetTitle("Boxes");
Int_t i,j;
Int_t nx = 40;
Int_t ny = 40;
Double_t xval1,yval1,xval2,yval2;
Double_t dx=0.2, dy=0.1;
xval1 = 0.;
xval2 = dx;
for (i = 0; i<nx; i++) {
yval1 = 0.;
yval2 = dy;
for (j = 0; j<ny; j++) {
h2p->AddBin(xval1, yval1, xval2, yval2);
yval1 = yval2;
yval2 = yval2+yval2*dy;
}
xval1 = xval2;
xval2 = xval2+xval2*dx;
}
TRandom ran;
for (i = 0; i<300000; i++) {
h2p->Fill(50*ran.Gaus(2.,1), ran.Gaus(2.,1));
}
h2p->Draw("COLZ");
return ch2p2;
}
virtual void SetTitle(const char *title)
See GetStatOverflows for more information.
Definition: TH1.cxx:6333
virtual void SetName(const char *name)
Change the name of this histogram.
Definition: TH1.cxx:8404
virtual Int_t Fill(Double_t x, Double_t y)
Increment the bin containing (x,y) by 1.
Definition: TH2Poly.cxx:589
virtual Int_t AddBin(TObject *poly)
Adds a new bin to the histogram.
Definition: TH2Poly.cxx:222

One TH2Poly bin can be a list of polygons. Such bins are defined by calling AddBin with a TMultiGraph. The following example shows a such case:

{
auto ch2p2 = new TCanvas("ch2p2","ch2p2",600,400);
Int_t i, bin;
const Int_t nx = 48;
const char *states [nx] = {
"alabama", "arizona", "arkansas", "california",
"colorado", "connecticut", "delaware", "florida",
"georgia", "idaho", "illinois", "indiana",
"iowa", "kansas", "kentucky", "louisiana",
"maine", "maryland", "massachusetts", "michigan",
"minnesota", "mississippi", "missouri", "montana",
"nebraska", "nevada", "new_hampshire", "new_jersey",
"new_mexico", "new_york", "north_carolina", "north_dakota",
"ohio", "oklahoma", "oregon", "pennsylvania",
"rhode_island", "south_carolina", "south_dakota", "tennessee",
"texas", "utah", "vermont", "virginia",
"washington", "west_virginia", "wisconsin", "wyoming"
};
Double_t pop[nx] = {
4708708, 6595778, 2889450, 36961664, 5024748, 3518288, 885122, 18537969,
9829211, 1545801, 12910409, 6423113, 3007856, 2818747, 4314113, 4492076,
1318301, 5699478, 6593587, 9969727, 5266214, 2951996, 5987580, 974989,
1796619, 2643085, 1324575, 8707739, 2009671, 19541453, 9380884, 646844,
11542645, 3687050, 3825657, 12604767, 1053209, 4561242, 812383, 6296254,
24782302, 2784572, 621760, 7882590, 6664195, 1819777, 5654774, 544270
};
Double_t lon1 = -130;
Double_t lon2 = -65;
Double_t lat1 = 24;
Double_t lat2 = 50;
auto p = new TH2Poly("USA","USA Population",lon1,lon2,lat1,lat2);
auto f = TFile::Open("http://root.cern.ch/files/usa.root", "CACHEREAD");
TKey *key;
TIter nextkey(gDirectory->GetListOfKeys());
while ((key = (TKey*)nextkey())) {
TObject *obj = key->ReadObj();
if (obj->InheritsFrom("TMultiGraph")) {
mg = (TMultiGraph*)obj;
bin = p->AddBin(mg);
}
}
for (i=0; i<nx; i++) p->Fill(states[i], pop[i]);
p->Draw("COLZ L");
}
static Bool_t SetCacheFileDir(ROOT::Internal::TStringView cacheDir, Bool_t operateDisconnected=kTRUE, Bool_t forceCacheread=kFALSE)
Definition: TFile.h:319
Book space in a file, create I/O buffers, to fill them, (un)compress them.
Definition: TKey.h:24
virtual TObject * ReadObj()
To read a TObject* from the file.
Definition: TKey.cxx:729
A TMultiGraph is a collection of TGraph (or derived) objects.
Definition: TMultiGraph.h:35
Mother of all ROOT objects.
Definition: TObject.h:37
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition: TObject.cxx:443
static constexpr double mg

TH2Poly histograms can also be plotted using the GL interface using the option "GLLEGO".

Since
ROOT version 6.09/01

In some cases it can be useful to not draw the empty bins. the option "0" combined with the option "COL" et COLZ allows to do that.

{
auto chc = new TCanvas("chc","chc",600,400);
auto hc = new TH2Poly();
hc->Honeycomb(0,0,.1,25,25);
hc->SetName("hc");
hc->SetTitle("Option COLZ 0");
TRandom ran;
for (int i = 0; i<300; i++) hc->Fill(ran.Gaus(2.,1), ran.Gaus(2.,1));
hc->Draw("colz 0");
}

The SPEC option

This option allows to use the TSpectrum2Painter tools. See the full documentation in TSpectrum2Painter::PaintSpectrum.

Option "Z" : Adding the color palette on the right side of the pad

When this option is specified, a color palette with an axis indicating the value of the corresponding color is drawn on the right side of the picture. In case, not enough space is left, one can increase the size of the right margin by calling TPad::SetRightMargin(). The attributes used to display the palette axis values are taken from the Z axis of the object. For example, to set the labels size on the palette axis do:

hist->GetZaxis()->SetLabelSize().

WARNING: The palette axis is always drawn vertically.

Setting the color palette

To change the color palette TStyle::SetPalette should be used, eg:

gStyle->SetPalette(ncolors,colors);

For example the option COL draws a 2D histogram with cells represented by a box filled with a color index which is a function of the cell content. If the cell content is N, the color index used will be the color number in colors[N], etc. If the maximum cell content is greater than ncolors, all cell contents are scaled to ncolors.

If ncolors <= 0, a default palette (see below) of 50 colors is defined. This palette is recommended for pads, labels ...

if ncolors == 1 && colors == 0, then a Pretty Palette with a Spectrum Violet->Red is created with 50 colors. That's the default rain bow palette.

Other pre-defined palettes with 255 colors are available when colors == 0. The following value of ncolors give access to:

if ncolors = 51 and colors=0, a Deep Sea palette is used.
if ncolors = 52 and colors=0, a Grey Scale palette is used.
if ncolors = 53 and colors=0, a Dark Body Radiator palette is used.
if ncolors = 54 and colors=0, a two-color hue palette palette is used.(dark blue through neutral gray to bright yellow)
if ncolors = 55 and colors=0, a Rain Bow palette is used.
if ncolors = 56 and colors=0, an inverted Dark Body Radiator palette is used.

If ncolors > 0 && colors == 0, the default palette is used with a maximum of ncolors.

The default palette defines:

  • index 0 to 9 : shades of grey
  • index 10 to 19 : shades of brown
  • index 20 to 29 : shades of blue
  • index 30 to 39 : shades of red
  • index 40 to 49 : basic colors

The color numbers specified in the palette can be viewed by selecting the item colors in the VIEW menu of the canvas tool bar. The red, green, and blue components of a color can be changed thanks to TColor::SetRGB().

Since
ROOT version 6.19/01

As default labels and ticks are drawn by TGAxis at equidistant (lin or log) points as controlled by SetNdivisions. If option "CJUST" is given labels and ticks are justified at the color boundaries defined by the contour levels. For more details see TPaletteAxis

Drawing a sub-range of a 2D histogram; the [cutg] option

Using a TCutG object, it is possible to draw a sub-range of a 2D histogram. One must create a graphical cut (mouse or C++) and specify the name of the cut between [] in the Draw() option. For example (fit2a.C), with a TCutG named cutg, one can call:

myhist->Draw("surf1 [cutg]");

To invert the cut, it is enough to put a - in front of its name:

myhist->Draw("surf1 [-cutg]");

It is possible to apply several cuts (, means logical AND):

myhist->Draw("surf1 [cutg1,cutg2]");
#include "TF2.h"
#include "TH2.h"
#include "TCutG.h"
#include "TMath.h"
#include "TCanvas.h"
#include "TStyle.h"
Double_t r1 = Double_t((x[0]-par[1])/par[2]);
Double_t r2 = Double_t((x[1]-par[3])/par[4]);
return par[0]*TMath::Exp(-0.5*(r1*r1+r2*r2));
}
Double_t fun2(Double_t *x, Double_t *par) {
Double_t *p1 = &par[0];
Double_t *p2 = &par[5];
Double_t *p3 = &par[10];
Double_t result = g2(x,p1) + g2(x,p2) + g2(x,p3);
return result;
}
TCanvas *fit2a() {
TCanvas *c = new TCanvas();
const Int_t npar = 15;
Double_t f2params[npar] = {100,-3,3,-3,3,160,0,0.8,0,0.9,40,4,0.7,4,0.7};
auto f2 = new TF2("f2",fun2,-10,10,-10,10, npar);
f2->SetParameters(f2params);
//Create an histogram and fill it randomly with f2
auto h2 = new TH2F("h2","From f2",40,-10,10,40,-10,10);
Int_t nentries = 100000;
h2->FillRandom("f2",nentries);
//Fit h2 with original function f2
Float_t ratio = 4*nentries/100000;
f2params[ 0] *= ratio;
f2params[ 5] *= ratio;
f2params[10] *= ratio;
f2->SetParameters(f2params);
h2->Fit("f2","N");
auto cutg = new TCutG("cutg",5);
cutg->SetPoint(0,-7,-7);
cutg->SetPoint(1, 2,-7);
cutg->SetPoint(2, 2, 2);
cutg->SetPoint(3,-7, 2);
cutg->SetPoint(4,-7,-7);
h2->Draw("lego2 0");
h2->SetFillColor(38);
f2->SetNpx(80);
f2->SetNpy(80);
f2->Draw("surf1 same bb [cutg]");
return c;
}
const Bool_t kTRUE
Definition: RtypesCore.h:87
int nentries
Definition: THbookFile.cxx:89
Graphical cut class.
Definition: TCutG.h:20
A 2-Dim function with parameters.
Definition: TF2.h:29
virtual void FillRandom(const char *fname, Int_t ntimes=5000)
Fill histogram following distribution in function fname.
Definition: TH1.cxx:3445
Double_t Exp(Double_t x)
Definition: TMath.h:717

Drawing options for 3D histograms

Option Description
"ISO" Draw a Gouraud shaded 3d iso surface through a 3d histogram. It paints one surface at the value computed as follow: SumOfWeights/(NbinsX*NbinsY*NbinsZ)
"BOX" Draw a for each cell with volume proportional to the content's absolute value. An hidden line removal algorithm is used
"BOX1" Same as BOX but an hidden surface removal algorithm is used
"BOX2" The boxes' colors are picked in the current palette according to the bins' contents
"BOX2Z" Same as "BOX2". In addition the color palette is also drawn.
"BOX3" Same as BOX1, but the border lines of each lego-bar are not drawn.

Note that instead of BOX one can also use LEGO.

By default, like 2D histograms, 3D histograms are drawn as scatter plots.

The following example shows a 3D histogram plotted as a scatter plot.

{
auto c06 = new TCanvas("c06","c06",600,400);
auto h3scat = new TH3F("h3scat","Option SCAT (default) ",15,-2,2,15,-2,2,15,0,4);
double x, y, z;
for (Int_t i=0;i<10000;i++) {
z = x*x + y*y;
h3scat->Fill(x,y,z);
}
h3scat->Draw();
}
3-D histogram with a float per channel (see TH1 documentation)}
Definition: TH3.h:267

The following example shows a 3D histogram plotted with the option BOX.

{
auto c16 = new TCanvas("c16","c16",600,400);
auto h3box = new TH3F("h3box","Option BOX",15,-2,2,15,-2,2,15,0,4);
double x, y, z;
for (Int_t i=0;i<10000;i++) {
z = x*x + y*y;
h3box->Fill(x,y,z);
}
h3box->Draw("BOX");
}

The following example shows a 3D histogram plotted with the option BOX1.

{
auto c36 = new TCanvas("c36","c36",600,400);
auto h3box = new TH3F("h3box","Option BOX1",10,-2.,2.,10,-2.,2.,10,-0.5,2.);
double x, y, z;
for (Int_t i=0;i<10000;i++) {
z = abs(sin(x)/x + cos(y)*y);
h3box->Fill(x,y,z);
}
h3box->SetFillColor(9);
h3box->Draw("BOX1");
}
double cos(double)

The following example shows a 3D histogram plotted with the option BOX2.

{
auto c56 = new TCanvas("c56","c56",600,400);
auto h3box = new TH3F("h3box","Option BOX2",10,-2.,2.,10,-2.,2.,10,-0.5,2.);
double x, y, z;
for (Int_t i=0;i<10000;i++) {
z = abs(sin(x)/x + cos(y)*y);
h3box->Fill(x,y,z);
}
h3box->Draw("BOX2 Z");
}

The following example shows a 3D histogram plotted with the option BOX3.

{
auto c46 = new TCanvas("c46","c46",600,400);
c46->SetFillColor(38);
auto h3box = new TH3F("h3box","Option BOX3",15,-2,2,15,-2,2,15,0,4);
double x, y, z;
for (Int_t i=0;i<10000;i++) {
z = x*x + y*y;
h3box->Fill(x,y,z);
}
h3box->Draw("BOX3");
}

For all the BOX options each bin is drawn as a 3D box with a volume proportional to the absolute value of the bin content. The bins with a negative content are drawn with a X on each face of the box as shown in the following example:

{
auto c = new TCanvas("c","c",600,400);
auto h3box = new TH3F("h3box","Option BOX1 with negative bins",3, 0., 4., 3, 0.,4., 3, 0., 4.);
h3box->Fill(0., 2., 2., 10.);
h3box->Fill(2., 2., 2., 5.);
h3box->Fill(2., 2., .5, 2.);
h3box->Fill(2., 2., 3., -1.);
h3box->Fill(3., 2., 2., -10.);
h3box->SetFillColor(8);
h3box->Draw("box1");
}

The following example shows a 3D histogram plotted with the option ISO.

{
auto c26 = new TCanvas("c26","c26",600,400);
auto h3iso = new TH3F("h3iso","Option ISO",15,-2,2,15,-2,2,15,0,4);
double x, y, z;
for (Int_t i=0;i<10000;i++) {
z = x*x + y*y;
h3iso->Fill(x,y,z);
}
h3iso->SetFillColor(kCyan);
h3iso->Draw("ISO");
}
@ kCyan
Definition: Rtypes.h:64

Drawing option for histograms' stacks

Stacks of histograms are managed with the THStack. A THStack is a collection of TH1 (or derived) objects. For painting only the THStack containing TH1 only or THStack containing TH2 only will be considered.

By default, histograms are shown stacked:

  1. The first histogram is paint.
  2. The the sum of the first and second, etc...

If the option NOSTACK is specified, the histograms are all paint in the same pad as if the option SAME had been specified. This allows to compute X and Y scales common to all the histograms, like TMultiGraph does for graphs.

If the option PADS is specified, the current pad/canvas is subdivided into a number of pads equal to the number of histograms and each histogram is paint into a separate pad.

The following example shows various types of stacks (hstack.C).

TCanvas *hstack() {
THStack *hs = new THStack("hs","Stacked 1D histograms");
//create three 1-d histograms
TH1F *h1st = new TH1F("h1st","test hstack",100,-4,4);
h1st->FillRandom("gaus",20000);
h1st->SetMarkerStyle(21);
hs->Add(h1st);
TH1F *h2st = new TH1F("h2st","test hstack",100,-4,4);
h2st->FillRandom("gaus",15000);
h2st->SetMarkerStyle(21);
hs->Add(h2st);
TH1F *h3st = new TH1F("h3st","test hstack",100,-4,4);
h3st->FillRandom("gaus",10000);
h3st->SetMarkerStyle(21);
hs->Add(h3st);
TCanvas *cst = new TCanvas("cst","stacked hists",10,10,700,700);
cst->Divide(2,2);
// in top left pad, draw the stack with defaults
cst->cd(1);
hs->Draw();
// in top right pad, draw the stack in non-stack mode
// and errors option
cst->cd(2);
gPad->SetGrid();
hs->Draw("nostack,e1p");
//in bottom left, draw in stack mode with "lego1" option
cst->cd(3);
gPad->SetFrameFillColor(17);
gPad->SetTheta(3.77);
gPad->SetPhi(2.9);
hs->Draw("lego1");
cst->cd(4);
//create two 2-D histograms and draw them in stack mode
gPad->SetFrameFillColor(17);
THStack *a = new THStack("a","Stacked 2D histograms");
TF2 *f1 = new TF2("f1",
"xygaus + xygaus(5) + xylandau(10)",-4,4,-4,4);
Double_t params1[] = {130,-1.4,1.8,1.5,1, 150,2,0.5,-2,0.5,
3600,-2,0.7,-3,0.3};
f1->SetParameters(params1);
TH2F *h2sta = new TH2F("h2sta","h2sta",20,-4,4,20,-4,4);
h2sta->SetFillColor(38);
h2sta->FillRandom("f1",4000);
TF2 *f2 = new TF2("f2","xygaus + xygaus(5)",-4,4,-4,4);
Double_t params2[] = {100,-1.4,1.9,1.1,2, 80,2,0.7,-2,0.5};
f2->SetParameters(params2);
TH2F *h2stb = new TH2F("h2stb","h2stb",20,-4,4,20,-4,4);
h2stb->SetFillColor(46);
h2stb->FillRandom("f2",3000);
a->Add(h2sta);
a->Add(h2stb);
a->Draw();
return cst;
}
virtual void SetMarkerColor(Color_t mcolor=1)
Set the marker color.
Definition: TAttMarker.h:38
TVirtualPad * cd(Int_t subpadnumber=0)
Set current canvas & pad.
Definition: TCanvas.cxx:696
virtual void FillRandom(const char *fname, Int_t ntimes=5000)
Fill histogram following distribution in function fname.
Definition: TH2.cxx:597
virtual void Draw(Option_t *chopt="")
Draw this multihist with its current attributes.
Definition: THStack.cxx:445
virtual void Add(TH1 *h, Option_t *option="")
add a new histogram to the list Only 1-d and 2-d histograms currently supported.
Definition: THStack.cxx:359
virtual void Divide(Int_t nx=1, Int_t ny=1, Float_t xmargin=0.01, Float_t ymargin=0.01, Int_t color=0)
Automatic pad generation by division.
Definition: TPad.cxx:1163
auto * a
Definition: textangle.C:12

The option nostackb allows to draw the histograms next to each other as bar charts:

{
auto cst0 = new TCanvas("cst0","cst0",600,400);
auto hs = new THStack("hs","Stacked 1D histograms: option #font[82]{\"nostackb\"}");
auto h1 = new TH1F("h1","h1",10,-4,4);
h1->FillRandom("gaus",20000);
hs->Add(h1);
auto h2 = new TH1F("h2","h2",10,-4,4);
h2->FillRandom("gaus",15000);
hs->Add(h2);
auto h3 = new TH1F("h3","h3",10,-4,4);
h3->FillRandom("gaus",10000);
h3->SetFillColor(kGreen);
hs->Add(h3);
hs->Draw("nostackb");
hs->GetXaxis()->SetNdivisions(-10);
cst0->SetGridx();
}
virtual void SetNdivisions(Int_t n=510, Bool_t optim=kTRUE)
Set the number of divisions for this axis.
Definition: TAttAxis.cxx:229
TAxis * GetXaxis() const
Get x axis of the histogram used to draw the stack.
Definition: THStack.cxx:601

If at least one of the histograms in the stack has errors, the whole stack is visualized by default with error bars. To visualize it without errors the option HIST should be used.

{
auto cst1 = new TCanvas("cst1","cst1",700,400);
cst1->Divide(2,1);
auto hst11 = new TH1F("hst11", "", 20, -10, 10);
hst11->Sumw2();
hst11->FillRandom("gaus", 1000);
hst11->SetFillColor(kViolet);
hst11->SetLineColor(kViolet);
auto hst12 = new TH1F("hst12", "", 20, -10, 10);
hst12->FillRandom("gaus", 500);
hst12->SetFillColor(kBlue);
hst12->SetLineColor(kBlue);
THStack st1("st1", "st1");
st1.Add(hst11);
st1.Add(hst12);
cst1->cd(1); st1.Draw();
cst1->cd(2); st1.Draw("hist");
}
@ kViolet
Definition: Rtypes.h:65

Drawing of 3D implicit functions

3D implicit functions (TF3) can be drawn as iso-surfaces. The implicit function f(x,y,z) = 0 is drawn in cartesian coordinates. In the following example the options "FB" and "BB" suppress the "Front Box" and "Back Box" around the plot.

{
auto c2 = new TCanvas("c2","c2",600,400);
auto f3 = new TF3("f3","sin(x*x+y*y+z*z-36)",-2,2,-2,2,-2,2);
f3->SetClippingBoxOn(0,0,0);
f3->SetFillColor(30);
f3->SetLineColor(15);
f3->Draw("FBBB");
}
A 3-Dim function with parameters.
Definition: TF3.h:28

Associated functions drawing

An associated function is created by TH1::Fit. More than on fitted function can be associated with one histogram (see TH1::Fit).

A TF1 object f1 can be added to the list of associated functions of an histogram h without calling TH1::Fit simply doing:

h->GetListOfFunctions()->Add(f1);

or

h->GetListOfFunctions()->Add(f1,someoption);

To retrieve a function by name from this list, do:

TF1 *f1 = (TF1*)h->GetListOfFunctions()->FindObject(name);

or

TF1 *f1 = h->GetFunction(name);

Associated functions are automatically painted when an histogram is drawn. To avoid the painting of the associated functions the option HIST should be added to the list of the options used to paint the histogram.

Drawing using OpenGL

The class TGLHistPainter allows to paint data set using the OpenGL 3D graphics library. The plotting options start with GL keyword. In addition, in order to inform canvases that OpenGL should be used to render 3D representations, the following option should be set:

gStyle->SetCanvasPreferGL(true);

General information: plot types and supported options

The following types of plots are provided:

For lego plots the supported options are:

Option Description
"GLLEGO" Draw a lego plot. It works also for TH2Poly.
"GLLEGO2" Bins with color levels.
"GLLEGO3" Cylindrical bars.

Lego painter in cartesian supports logarithmic scales for X, Y, Z. In polar only Z axis can be logarithmic, in cylindrical only Y.

For surface plots (TF2 and TH2) the supported options are:

Option Description
"GLSURF" Draw a surface.
"GLSURF1" Surface with color levels
"GLSURF2" The same as "GLSURF1" but without polygon outlines.
"GLSURF3" Color level projection on top of plot (works only in cartesian coordinate system).
"GLSURF4" Same as "GLSURF" but without polygon outlines.

The surface painting in cartesian coordinates supports logarithmic scales along X, Y, Z axis. In polar coordinates only the Z axis can be logarithmic, in cylindrical coordinates only the Y axis.

Additional options to SURF and LEGO - Coordinate systems:

Option Description
" " Default, cartesian coordinates system.
"POL" Polar coordinates system.
"CYL" Cylindrical coordinates system.
"SPH" Spherical coordinates system.

TH3 as color boxes

The supported option is:

Option Description
"GLCOL" H3 is drawn using semi-transparent colored boxes. See $ROOTSYS/tutorials/gl/glvox1.C.

TH3 as boxes (spheres)

The supported options are:

Option Description
"GLBOX" TH3 as a set of boxes, size of box is proportional to bin content.
"GLBOX1" The same as "glbox", but spheres are drawn instead of boxes.

TH3 as iso-surface(s)

The supported option is:

Option Description
"GLISO" TH3 is drawn using iso-surfaces.

TF3 (implicit function)

The supported option is:

Option Description
"GLTF3" Draw a TF3.

Parametric surfaces

$ROOTSYS/tutorials/gl/glparametric.C shows how to create parametric equations and visualize the surface.

Interaction with the plots

All the interactions are implemented via standard methods DistancetoPrimitive() and ExecuteEvent(). That's why all the interactions with the OpenGL plots are possible only when the mouse cursor is in the plot's area (the plot's area is the part of a the pad occupied by gl-produced picture). If the mouse cursor is not above gl-picture, the standard pad interaction is performed.

Selectable parts

Different parts of the plot can be selected:

  • xoz, yoz, xoy back planes: When such a plane selected, it's highlighted in green if the dynamic slicing by this plane is supported, and it's highlighted in red, if the dynamic slicing is not supported.
  • The plot itself: On surfaces, the selected surface is outlined in red. (TF3 and ISO are not outlined). On lego plots, the selected bin is highlighted. The bin number and content are displayed in pad's status bar. In box plots, the box or sphere is highlighted and the bin info is displayed in pad's status bar.

Rotation and zooming

  • Rotation: When the plot is selected, it can be rotated by pressing and holding the left mouse button and move the cursor.
  • Zoom/Unzoom: Mouse wheel or 'j', 'J', 'k', 'K' keys.

Panning

The selected plot can be moved in a pad's area by pressing and holding the left mouse button and the shift key.

Box cut

Surface, iso, box, TF3 and parametric painters support box cut by pressing the 'c' or 'C' key when the mouse cursor is in a plot's area. That will display a transparent box, cutting away part of the surface (or boxes) in order to show internal part of plot. This box can be moved inside the plot's area (the full size of the box is equal to the plot's surrounding box) by selecting one of the box cut axes and pressing the left mouse button to move it.

Plot specific interactions (dynamic slicing etc.)

Currently, all gl-plots support some form of slicing. When back plane is selected (and if it's highlighted in green) you can press and hold left mouse button and shift key and move this back plane inside plot's area, creating the slice. During this "slicing" plot becomes semi-transparent. To remove all slices (and projected curves for surfaces) double click with left mouse button in a plot's area.

Surface with option "GLSURF"

The surface profile is displayed on the slicing plane. The profile projection is drawn on the back plane by pressing ‘'p’or'P'` key.

TF3

The contour plot is drawn on the slicing plane. For TF3 the color scheme can be changed by pressing 's' or 'S'.

Box

The contour plot corresponding to slice plane position is drawn in real time.

Iso

Slicing is similar to "GLBOX" option.

Parametric plot

No slicing. Additional keys: 's' or 'S' to change color scheme - about 20 color schemes supported ('s' for "scheme"); 'l' or 'L' to increase number of polygons ('l' for "level" of details), 'w' or 'W' to show outlines ('w' for "wireframe").

Highlight mode for histogram

Since
ROOT version 6.15/01
Highlight mode

Highlight mode is implemented for TH1 (and for TGraph) class. When highlight mode is on, mouse movement over the bin will be represented graphically. Bin will be highlighted as "bin box" (presented by box object). Moreover, any highlight (change of bin) emits signal TCanvas::Highlighted() which allows the user to react and call their own function. For a better understanding see also the tutorials $ROOTSYS/tutorials/hist/hlHisto*.C files.

Highlight mode is switched on/off by TH1::SetHighlight() function or interactively from TH1 context menu. TH1::IsHighlight() to verify whether the highlight mode enabled or disabled, default it is disabled.

root [0] .x $ROOTSYS/tutorials/hsimple.C
root [1] hpx->SetHighlight(kTRUE) // or interactively from TH1 context menu
root [2] hpx->IsHighlight()
(bool) true
Highlight mode for histogram

Highlight mode and user function

The user can use (connect) TCanvas::Highlighted() signal, which is always emitted if there is a highlight bin and call user function via signal and slot communication mechanism. TCanvas::Highlighted() is similar TCanvas::Picked()

  • when selected object (histogram as a whole) is different from previous then emit Picked() signal
  • when selected (highlighted) bin from histogram is different from previous then emit Highlighted() signal

Any user function (or functions) has to be defined UserFunction(TVirtualPad *pad, TObject *obj, Int_t x, Int_t y). In example (see below) has name PrintInfo(). All parameters of user function are taken from

void TCanvas::Highlighted(TVirtualPad *pad, TObject *obj, Int_t x, Int_t y)
  • pad is pointer to pad with highlighted histogram
  • obj is pointer to highlighted histogram
  • x is highlighted x bin for 1D histogram
  • y is highlighted y bin for 2D histogram (for 1D histogram not in use)

Example how to create a connection from any TCanvas object to a user UserFunction() slot (see also TQObject::Connect() for additional info)

TQObject::Connect("TCanvas", "Highlighted(TVirtualPad*,TObject*,Int_t,Int_t)",
                      0, 0, "UserFunction(TVirtualPad*,TObject*,Int_t,Int_t)");

or use non-static "simplified" function TCanvas::HighlightConnect(const char *slot)

c1->HighlightConnect("UserFunction(TVirtualPad*,TObject*,Int_t,Int_t)");

NOTE the signal and slot string must have a form "(TVirtualPad*,TObject*,Int_t,Int_t)"

root [0] .x $ROOTSYS/tutorials/hsimple.C
root [1] hpx->SetHighlight(kTRUE)
root [2] .x hlprint.C

file hlprint.C

void PrintInfo(TVirtualPad *pad, TObject *obj, Int_t x, Int_t y)
{
auto h = (TH1F *)obj;
if (!h->IsHighlight()) // after highlight disabled
h->SetTitle("highlight disable");
else
h->SetTitle(TString::Format("bin[%03d] (%5.2f) content %g", x,
h->GetBinCenter(x), h->GetBinContent(x)));
pad->Update();
}
void hlprint()
{
if (!gPad) return;
gPad->GetCanvas()->HighlightConnect("PrintInfo(TVirtualPad*,TObject*,Int_t,Int_t)");
}
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:2311
TVirtualPad is an abstract base class for the Pad and Canvas classes.
Definition: TVirtualPad.h:50
virtual void Update()=0
Highlight mode and simple user function

For more complex demo please see for example $ROOTSYS/tutorials/tree/temperature.C file.

Definition at line 47 of file THistPainter.h.

Public Member Functions

 THistPainter ()
 Default constructor. More...
 
virtual ~THistPainter ()
 Default destructor. More...
 
virtual std::vector< THistRenderingRegionComputeRenderingRegions (TAxis *pAxis, Int_t nPixels, bool isLog)
 Returns the rendering regions for an axis to use in the COL2 option. More...
 
virtual void DefineColorLevels (Int_t ndivz)
 Define the color levels used to paint legos, surfaces etc.. More...
 
virtual Int_t DistancetoPrimitive (Int_t px, Int_t py)
 Compute the distance from the point px,py to a line. More...
 
virtual void DrawPanel ()
 Display a panel with all histogram drawing options. More...
 
virtual void ExecuteEvent (Int_t event, Int_t px, Int_t py)
 Execute the actions corresponding to event. More...
 
virtual TListGetContourList (Double_t contour) const
 Get a contour (as a list of TGraphs) using the Delaunay triangulation. More...
 
virtual char * GetObjectInfo (Int_t px, Int_t py) const
 Display the histogram info (bin number, contents, integral up to bin corresponding to cursor position px,py. More...
 
virtual TListGetStack () const
 
virtual Int_t GetXHighlightBin () const
 
virtual Int_t GetYHighlightBin () const
 
virtual void HighlightBin (Int_t px, Int_t py)
 Check on highlight bin. More...
 
virtual Bool_t IsInside (Double_t x, Double_t y)
 Return kTRUE if the point x, y is inside one of the graphical cuts. More...
 
virtual Bool_t IsInside (Int_t x, Int_t y)
 Return kTRUE if the cell ix, iy is inside one of the graphical cuts. More...
 
virtual Int_t MakeChopt (Option_t *option)
 Decode string choptin and fill Hoption structure. More...
 
virtual Int_t MakeCuts (char *cutsopt)
 Decode string choptin and fill Graphical cuts structure. More...
 
virtual void Paint (Option_t *option="")
 Control routine to paint any kind of histograms More...
 
virtual void Paint2DErrors (Option_t *option)
 Draw 2D histograms errors. More...
 
virtual void PaintArrows (Option_t *option)
 Control function to draw a table as an arrow plot More...
 
virtual void PaintAxis (Bool_t drawGridOnly=kFALSE)
 Draw axis (2D case) of an histogram. More...
 
virtual void PaintBar (Option_t *option)
 Draw a bar-chart in a normal pad. More...
 
virtual void PaintBarH (Option_t *option)
 Draw a bar char in a rotated pad (X vertical, Y horizontal) More...
 
virtual void PaintBoxes (Option_t *option)
 Control function to draw a 2D histogram as a box plot More...
 
virtual void PaintCandlePlot (Option_t *option)
 Control function to draw a 2D histogram as a candle (box) plot or violin plot More...
 
virtual void PaintColorLevels (Option_t *option)
 Control function to draw a 2D histogram as a color plot. More...
 
virtual void PaintColorLevelsFast (Option_t *option)
 [Rendering scheme for the COL2 and COLZ2 options] (#HP14) More...
 
virtual void PaintContour (Option_t *option)
 Control function to draw a 2D histogram as a contour plot. More...
 
virtual Int_t PaintContourLine (Double_t elev1, Int_t icont1, Double_t x1, Double_t y1, Double_t elev2, Int_t icont2, Double_t x2, Double_t y2, Double_t *xarr, Double_t *yarr, Int_t *itarr, Double_t *levels)
 Fill the matrix xarr and yarr for Contour Plot. More...
 
virtual void PaintErrors (Option_t *option)
 Draw 1D histograms error bars. More...
 
virtual void PaintFrame ()
 Calculate range and clear pad (canvas). More...
 
virtual void PaintFunction (Option_t *option)
 [Paint functions associated to an histogram.](#HP28") More...
 
virtual void PaintH3 (Option_t *option="")
 Control function to draw a 3D histograms. More...
 
virtual void PaintH3Box (Int_t iopt)
 Control function to draw a 3D histogram with boxes. More...
 
virtual void PaintH3BoxRaster ()
 Control function to draw a 3D histogram with boxes. More...
 
virtual void PaintH3Iso ()
 Control function to draw a 3D histogram with Iso Surfaces. More...
 
virtual void PaintHighlightBin (Option_t *option="")
 Paint highlight bin as TBox object. More...
 
virtual void PaintHist (Option_t *option)
 Control routine to draw 1D histograms More...
 
virtual Int_t PaintInit ()
 Compute histogram parameters used by the drawing routines. More...
 
virtual Int_t PaintInitH ()
 Compute histogram parameters used by the drawing routines for a rotated pad. More...
 
virtual void PaintLego (Option_t *option)
 Control function to draw a 2D histogram as a lego plot. More...
 
virtual void PaintLegoAxis (TGaxis *axis, Double_t ang)
 Draw the axis for legos and surface plots. More...
 
virtual void PaintPalette ()
 Paint the color palette on the right side of the pad. More...
 
virtual void PaintScatterPlot (Option_t *option)
 Control function to draw a 2D histogram as a scatter plot. More...
 
virtual void PaintStat (Int_t dostat, TF1 *fit)
 Draw the statistics box for 1D and profile histograms. More...
 
virtual void PaintStat2 (Int_t dostat, TF1 *fit)
 Draw the statistics box for 2D histograms. More...
 
virtual void PaintStat3 (Int_t dostat, TF1 *fit)
 Draw the statistics box for 3D histograms. More...
 
virtual void PaintSurface (Option_t *option)
 Control function to draw a 2D histogram as a surface plot. More...
 
virtual void PaintTable (Option_t *option)
 Control function to draw 2D/3D histograms (tables). More...
 
virtual void PaintText (Option_t *option)
 Control function to draw a 1D/2D histograms with the bin values. More...
 
virtual void PaintTF3 ()
 Control function to draw a 3D implicit functions. More...
 
virtual void PaintTH2PolyBins (Option_t *option)
 Control function to draw a TH2Poly bins' contours. More...
 
virtual void PaintTH2PolyColorLevels (Option_t *option)
 Control function to draw a TH2Poly as a color plot. More...
 
virtual void PaintTH2PolyScatterPlot (Option_t *option)
 Control function to draw a TH2Poly as a scatter plot. More...
 
virtual void PaintTH2PolyText (Option_t *option)
 Control function to draw a TH2Poly as a text plot. More...
 
virtual void PaintTitle ()
 Draw the histogram title. More...
 
virtual void PaintTriangles (Option_t *option)
 Control function to draw a table using Delaunay triangles. More...
 
virtual void ProcessMessage (const char *mess, const TObject *obj)
 Process message mess. More...
 
virtual void RecalculateRange ()
 Recompute the histogram range following graphics operations. More...
 
virtual void RecursiveRemove (TObject *)
 Recursively remove this object from a list. More...
 
virtual void SetHighlight ()
 Set highlight (enable/disable) mode for fH. More...
 
virtual void SetHistogram (TH1 *h)
 Set current histogram to h More...
 
virtual void SetShowProjection (const char *option, Int_t nbins)
 Set projection. More...
 
virtual void SetStack (TList *stack)
 
virtual void ShowProjection3 (Int_t px, Int_t py)
 Show projection (specified by fShowProjection) of a TH3. More...
 
virtual void ShowProjectionX (Int_t px, Int_t py)
 Show projection onto X. More...
 
virtual void ShowProjectionY (Int_t px, Int_t py)
 Show projection onto Y. More...
 
virtual Int_t TableInit ()
 Initialize various options to draw 2D histograms. More...
 
- Public Member Functions inherited from TVirtualHistPainter
 TVirtualHistPainter ()
 
virtual ~TVirtualHistPainter ()
 
virtual Int_t DistancetoPrimitive (Int_t px, Int_t py)=0
 Computes distance from point (px,py) to the object. More...
 
virtual void DrawPanel ()=0
 
virtual void ExecuteEvent (Int_t event, Int_t px, Int_t py)=0
 Execute action corresponding to an event at (px,py). More...
 
virtual TListGetContourList (Double_t contour) const =0
 
virtual char * GetObjectInfo (Int_t px, Int_t py) const =0
 Returns string containing info about the object at position (px,py). More...
 
virtual TListGetStack () const =0
 
virtual Bool_t IsInside (Double_t x, Double_t y)=0
 
virtual Bool_t IsInside (Int_t x, Int_t y)=0
 
virtual Int_t MakeCuts (char *cutsopt)=0
 
virtual void Paint (Option_t *option="")=0
 This method must be overridden if a class wants to paint itself. More...
 
virtual void PaintStat (Int_t dostat, TF1 *fit)=0
 
virtual void ProcessMessage (const char *mess, const TObject *obj)=0
 
virtual void SetHighlight ()=0
 
virtual void SetHistogram (TH1 *h)=0
 
virtual void SetShowProjection (const char *option, Int_t nbins)=0
 
virtual void SetStack (TList *stack)=0
 
- Public Member Functions inherited from TObject
 TObject ()
 TObject constructor. More...
 
 TObject (const TObject &object)
 TObject copy ctor. More...
 
virtual ~TObject ()
 TObject destructor. More...
 
void AbstractMethod (const char *method) const
 Use this method to implement an "abstract" method that you don't want to leave purely abstract. More...
 
virtual void AppendPad (Option_t *option="")
 Append graphics object to current pad. More...
 
virtual void Browse (TBrowser *b)
 Browse object. May be overridden for another default action. More...
 
ULong_t CheckedHash ()
 Check and record whether this class has a consistent Hash/RecursiveRemove setup (*) and then return the regular Hash value for this object. More...
 
virtual const char * ClassName () const
 Returns name of class to which the object belongs. More...
 
virtual void Clear (Option_t *="")
 
virtual TObjectClone (const char *newname="") const
 Make a clone of an object using the Streamer facility. More...
 
virtual Int_t Compare (const TObject *obj) const
 Compare abstract method. More...
 
virtual void Copy (TObject &object) const
 Copy this to obj. More...
 
virtual void Delete (Option_t *option="")
 Delete this object. More...
 
virtual Int_t DistancetoPrimitive (Int_t px, Int_t py)
 Computes distance from point (px,py) to the object. More...
 
virtual void Draw (Option_t *option="")
 Default Draw method for all objects. More...
 
virtual void DrawClass () const
 Draw class inheritance tree of the class to which this object belongs. More...
 
virtual TObjectDrawClone (Option_t *option="") const
 Draw a clone of this object in the current selected pad for instance with: gROOT->SetSelectedPad(gPad). More...
 
virtual void Dump () const
 Dump contents of object on stdout. More...
 
virtual void Error (const char *method, const char *msgfmt,...) const
 Issue error message. More...
 
virtual void Execute (const char *method, const char *params, Int_t *error=0)
 Execute method on this object with the given parameter string, e.g. More...
 
virtual void Execute (TMethod *method, TObjArray *params, Int_t *error=0)
 Execute method on this object with parameters stored in the TObjArray. More...
 
virtual void ExecuteEvent (Int_t event, Int_t px, Int_t py)
 Execute action corresponding to an event at (px,py). More...
 
virtual void Fatal (const char *method, const char *msgfmt,...) const
 Issue fatal error message. More...
 
virtual TObjectFindObject (const char *name) const
 Must be redefined in derived classes. More...
 
virtual TObjectFindObject (const TObject *obj) const
 Must be redefined in derived classes. More...
 
virtual Option_tGetDrawOption () const
 Get option used by the graphics system to draw this object. More...
 
virtual const char * GetIconName () const
 Returns mime type name of object. More...
 
virtual const char * GetName () const
 Returns name of object. More...
 
virtual char * GetObjectInfo (Int_t px, Int_t py) const
 Returns string containing info about the object at position (px,py). More...
 
virtual Option_tGetOption () const
 
virtual const char * GetTitle () const
 Returns title of object. More...
 
virtual UInt_t GetUniqueID () const
 Return the unique object id. More...
 
virtual Bool_t HandleTimer (TTimer *timer)
 Execute action in response of a timer timing out. More...
 
virtual ULong_t Hash () const
 Return hash value for this object. More...
 
Bool_t HasInconsistentHash () const
 Return true is the type of this object is known to have an inconsistent setup for Hash and RecursiveRemove (i.e. More...
 
virtual void Info (const char *method, const char *msgfmt,...) const
 Issue info message. More...
 
virtual Bool_t InheritsFrom (const char *classname) const
 Returns kTRUE if object inherits from class "classname". More...
 
virtual Bool_t InheritsFrom (const TClass *cl) const
 Returns kTRUE if object inherits from TClass cl. More...
 
virtual void Inspect () const
 Dump contents of this object in a graphics canvas. More...
 
void InvertBit (UInt_t f)
 
virtual Bool_t IsEqual (const TObject *obj) const
 Default equal comparison (objects are equal if they have the same address in memory). More...
 
virtual Bool_t IsFolder () const
 Returns kTRUE in case object contains browsable objects (like containers or lists of other objects). More...
 
R__ALWAYS_INLINE Bool_t IsOnHeap () const
 
virtual Bool_t IsSortable () const
 
R__ALWAYS_INLINE Bool_t IsZombie () const
 
virtual void ls (Option_t *option="") const
 The ls function lists the contents of a class on stdout. More...
 
void MayNotUse (const char *method) const
 Use this method to signal that a method (defined in a base class) may not be called in a derived class (in principle against good design since a child class should not provide less functionality than its parent, however, sometimes it is necessary). More...
 
virtual Bool_t Notify ()
 This method must be overridden to handle object notification. More...
 
void Obsolete (const char *method, const char *asOfVers, const char *removedFromVers) const
 Use this method to declare a method obsolete. More...
 
void operator delete (void *ptr)
 Operator delete. More...
 
void operator delete[] (void *ptr)
 Operator delete []. More...
 
voidoperator new (size_t sz)
 
voidoperator new (size_t sz, void *vp)
 
voidoperator new[] (size_t sz)
 
voidoperator new[] (size_t sz, void *vp)
 
TObjectoperator= (const TObject &rhs)
 TObject assignment operator. More...
 
virtual void Paint (Option_t *option="")
 This method must be overridden if a class wants to paint itself. More...
 
virtual void Pop ()
 Pop on object drawn in a pad to the top of the display list. More...
 
virtual void Print (Option_t *option="") const
 This method must be overridden when a class wants to print itself. More...
 
virtual Int_t Read (const char *name)
 Read contents of object with specified name from the current directory. More...
 
virtual void RecursiveRemove (TObject *obj)
 Recursively remove this object from a list. More...
 
void ResetBit (UInt_t f)
 
virtual void SaveAs (const char *filename="", Option_t *option="") const
 Save this object in the file specified by filename. More...
 
virtual void SavePrimitive (std::ostream &out, Option_t *option="")
 Save a primitive as a C++ statement(s) on output stream "out". More...
 
void SetBit (UInt_t f)
 
void SetBit (UInt_t f, Bool_t set)
 Set or unset the user status bits as specified in f. More...
 
virtual void SetDrawOption (Option_t *option="")
 Set drawing option for object. More...
 
virtual void SetUniqueID (UInt_t uid)
 Set the unique object id. More...
 
virtual void SysError (const char *method, const char *msgfmt,...) const
 Issue system error message. More...
 
R__ALWAYS_INLINE Bool_t TestBit (UInt_t f) const
 
Int_t TestBits (UInt_t f) const
 
virtual void UseCurrentStyle ()
 Set current style settings in this object This function is called when either TCanvas::UseCurrentStyle or TROOT::ForceStyle have been invoked. More...
 
virtual void Warning (const char *method, const char *msgfmt,...) const
 Issue warning message. More...
 
virtual Int_t Write (const char *name=0, Int_t option=0, Int_t bufsize=0)
 Write this object to the current directory. More...
 
virtual Int_t Write (const char *name=0, Int_t option=0, Int_t bufsize=0) const
 Write this object to the current directory. More...
 

Static Public Member Functions

static const char * GetBestFormat (Double_t v, Double_t e, const char *f)
 This function returns the best format to print the error value (e) knowing the parameter value (v) and the format (f) used to print it. More...
 
static void PaintSpecialObjects (const TObject *obj, Option_t *option)
 Static function to paint special objects like vectors and matrices. More...
 
static Int_t ProjectAitoff2xy (Double_t l, Double_t b, Double_t &Al, Double_t &Ab)
 Static function. More...
 
static Int_t ProjectMercator2xy (Double_t l, Double_t b, Double_t &Al, Double_t &Ab)
 Static function. More...
 
static Int_t ProjectParabolic2xy (Double_t l, Double_t b, Double_t &Al, Double_t &Ab)
 Static function code from Ernst-Jan Buis. More...
 
static Int_t ProjectSinusoidal2xy (Double_t l, Double_t b, Double_t &Al, Double_t &Ab)
 Static function code from Ernst-Jan Buis. More...
 
- Static Public Member Functions inherited from TVirtualHistPainter
static TVirtualHistPainterHistPainter (TH1 *obj)
 Static function returning a pointer to the current histogram painter. More...
 
static void SetPainter (const char *painter)
 Static function to set an alternative histogram painter. More...
 
- Static Public Member Functions inherited from TObject
static Long_t GetDtorOnly ()
 Return destructor only flag. More...
 
static Bool_t GetObjectStat ()
 Get status of object stat flag. More...
 
static void SetDtorOnly (void *obj)
 Set destructor only flag. More...
 
static void SetObjectStat (Bool_t stat)
 Turn on/off tracking of objects in the TObjectTable. More...
 

Protected Attributes

TCutGfCuts [kMaxCuts]
 
Int_t fCutsOpt [kMaxCuts]
 
TListfFunctions
 
TGraph2DPainterfGraph2DPainter
 
TH1fH
 
TPainter3dAlgorithmsfLego
 
Int_t fNcuts
 
TPiefPie
 
TString fShowOption
 
Int_t fShowProjection
 
TListfStack
 
TAxisfXaxis
 
Double_tfXbuf
 
Int_t fXHighlightBin
 
TAxisfYaxis
 
Double_tfYbuf
 
Int_t fYHighlightBin
 
TAxisfZaxis
 

Private Attributes

TString fObjectInfo
 

Additional Inherited Members

- Public Types inherited from TObject
enum  {
  kIsOnHeap = 0x01000000 , kNotDeleted = 0x02000000 , kZombie = 0x04000000 , kInconsistent = 0x08000000 ,
  kBitMask = 0x00ffffff
}
 
enum  { kSingleKey = BIT(0) , kOverwrite = BIT(1) , kWriteDelete = BIT(2) }
 
enum  EDeprecatedStatusBits { kObjInCanvas = BIT(3) }
 
enum  EStatusBits {
  kCanDelete = BIT(0) , kMustCleanup = BIT(3) , kIsReferenced = BIT(4) , kHasUUID = BIT(5) ,
  kCannotPick = BIT(6) , kNoContextMenu = BIT(8) , kInvalidObject = BIT(13)
}
 
- Protected Member Functions inherited from TObject
virtual void DoError (int level, const char *location, const char *fmt, va_list va) const
 Interface to ErrorHandler (protected). More...
 
void MakeZombie ()
 

#include <THistPainter.h>

Inheritance diagram for THistPainter:
[legend]

Constructor & Destructor Documentation

◆ THistPainter()

THistPainter::THistPainter ( )

Default constructor.

Definition at line 3123 of file THistPainter.cxx.

◆ ~THistPainter()

THistPainter::~THistPainter ( )
virtual

Default destructor.

Definition at line 3173 of file THistPainter.cxx.

Member Function Documentation

◆ ComputeRenderingRegions()

std::vector< THistRenderingRegion > THistPainter::ComputeRenderingRegions ( TAxis pAxis,
Int_t  nPixels,
bool  isLog 
)
virtual

Returns the rendering regions for an axis to use in the COL2 option.

The algorithm analyses the size of the axis compared to the size of the rendering region. It figures out the boundaries to use for each color of the rendering region. Only one axis is computed here.

This allows for a single computation of the boundaries before iterating through all of the bins.

Parameters
pAxisthe axis to consider
nPixelsthe number of pixels to render axis into
isLogwhether the axis is log scale

Definition at line 5373 of file THistPainter.cxx.

◆ DefineColorLevels()

void THistPainter::DefineColorLevels ( Int_t  ndivz)
virtual

Define the color levels used to paint legos, surfaces etc..

Definition at line 9447 of file THistPainter.cxx.

◆ DistancetoPrimitive()

Int_t THistPainter::DistancetoPrimitive ( Int_t  px,
Int_t  py 
)
virtual

Compute the distance from the point px,py to a line.

Compute the closest distance of approach from point px,py to elements of an histogram. The distance is computed in pixels units.

Algorithm: Currently, this simple model computes the distance from the mouse to the histogram contour only.

Implements TVirtualHistPainter.

Definition at line 3186 of file THistPainter.cxx.

◆ DrawPanel()

void THistPainter::DrawPanel ( )
virtual

Display a panel with all histogram drawing options.

Implements TVirtualHistPainter.

Definition at line 3374 of file THistPainter.cxx.

◆ ExecuteEvent()

void THistPainter::ExecuteEvent ( Int_t  event,
Int_t  px,
Int_t  py 
)
virtual

Execute the actions corresponding to event.

This function is called when a histogram is clicked with the locator at the pixel position px,py.

Implements TVirtualHistPainter.

Definition at line 3394 of file THistPainter.cxx.

◆ GetBestFormat()

const char * THistPainter::GetBestFormat ( Double_t  v,
Double_t  e,
const char *  f 
)
static

This function returns the best format to print the error value (e) knowing the parameter value (v) and the format (f) used to print it.

Definition at line 10561 of file THistPainter.cxx.

◆ GetContourList()

TList * THistPainter::GetContourList ( Double_t  contour) const
virtual

Get a contour (as a list of TGraphs) using the Delaunay triangulation.

Implements TVirtualHistPainter.

Definition at line 3634 of file THistPainter.cxx.

◆ GetObjectInfo()

char * THistPainter::GetObjectInfo ( Int_t  px,
Int_t  py 
) const
virtual

Display the histogram info (bin number, contents, integral up to bin corresponding to cursor position px,py.

Implements TVirtualHistPainter.

Definition at line 3662 of file THistPainter.cxx.

◆ GetStack()

virtual TList * THistPainter::GetStack ( ) const
inlinevirtual

Implements TVirtualHistPainter.

Definition at line 81 of file THistPainter.h.

◆ GetXHighlightBin()

virtual Int_t THistPainter::GetXHighlightBin ( ) const
inlinevirtual

Definition at line 82 of file THistPainter.h.

◆ GetYHighlightBin()

virtual Int_t THistPainter::GetYHighlightBin ( ) const
inlinevirtual

Definition at line 83 of file THistPainter.h.

◆ HighlightBin()

void THistPainter::HighlightBin ( Int_t  px,
Int_t  py 
)
virtual

Check on highlight bin.

Definition at line 3806 of file THistPainter.cxx.

◆ IsInside() [1/2]

Bool_t THistPainter::IsInside ( Double_t  x,
Double_t  y 
)
virtual

Return kTRUE if the point x, y is inside one of the graphical cuts.

Implements TVirtualHistPainter.

Definition at line 3945 of file THistPainter.cxx.

◆ IsInside() [2/2]

Bool_t THistPainter::IsInside ( Int_t  x,
Int_t  y 
)
virtual

Return kTRUE if the cell ix, iy is inside one of the graphical cuts.

Implements TVirtualHistPainter.

Definition at line 3927 of file THistPainter.cxx.

◆ MakeChopt()

Int_t THistPainter::MakeChopt ( Option_t option)
virtual

Decode string choptin and fill Hoption structure.

Definition at line 3961 of file THistPainter.cxx.

◆ MakeCuts()

Int_t THistPainter::MakeCuts ( char *  cutsopt)
virtual

Decode string choptin and fill Graphical cuts structure.

Implements TVirtualHistPainter.

Definition at line 4342 of file THistPainter.cxx.

◆ Paint()

void THistPainter::Paint ( Option_t option = "")
virtual

Control routine to paint any kind of histograms

Implements TVirtualHistPainter.

Definition at line 4389 of file THistPainter.cxx.

◆ Paint2DErrors()

void THistPainter::Paint2DErrors ( Option_t option)
virtual

Draw 2D histograms errors.

Definition at line 6521 of file THistPainter.cxx.

◆ PaintArrows()

void THistPainter::PaintArrows ( Option_t option)
virtual

Control function to draw a table as an arrow plot

Definition at line 4578 of file THistPainter.cxx.

◆ PaintAxis()

void THistPainter::PaintAxis ( Bool_t  drawGridOnly = kFALSE)
virtual

Draw axis (2D case) of an histogram.

If drawGridOnly is TRUE, only the grid is painted (if needed). This allows to draw the grid and the axis separately. In THistPainter::Paint this feature is used to make sure that the grid is drawn in the background and the axis tick marks in the foreground of the pad.

Definition at line 4676 of file THistPainter.cxx.

◆ PaintBar()

void THistPainter::PaintBar ( Option_t option)
virtual

Draw a bar-chart in a normal pad.

Definition at line 4953 of file THistPainter.cxx.

◆ PaintBarH()

void THistPainter::PaintBarH ( Option_t option)
virtual

Draw a bar char in a rotated pad (X vertical, Y horizontal)

Definition at line 5002 of file THistPainter.cxx.

◆ PaintBoxes()

void THistPainter::PaintBoxes ( Option_t option)
virtual

Control function to draw a 2D histogram as a box plot

Definition at line 5079 of file THistPainter.cxx.

◆ PaintCandlePlot()

void THistPainter::PaintCandlePlot ( Option_t option)
virtual

Control function to draw a 2D histogram as a candle (box) plot or violin plot

Definition at line 5264 of file THistPainter.cxx.

◆ PaintColorLevels()

void THistPainter::PaintColorLevels ( Option_t option)
virtual

Control function to draw a 2D histogram as a color plot.

Definition at line 5667 of file THistPainter.cxx.

◆ PaintColorLevelsFast()

void THistPainter::PaintColorLevelsFast ( Option_t option)
virtual

[Rendering scheme for the COL2 and COLZ2 options] (#HP14)

Definition at line 5478 of file THistPainter.cxx.

◆ PaintContour()

void THistPainter::PaintContour ( Option_t option)
virtual

Control function to draw a 2D histogram as a contour plot.

Definition at line 5830 of file THistPainter.cxx.

◆ PaintContourLine()

Int_t THistPainter::PaintContourLine ( Double_t  elev1,
Int_t  icont1,
Double_t  x1,
Double_t  y1,
Double_t  elev2,
Int_t  icont2,
Double_t  x2,
Double_t  y2,
Double_t xarr,
Double_t yarr,
Int_t itarr,
Double_t levels 
)
virtual

Fill the matrix xarr and yarr for Contour Plot.

Definition at line 6171 of file THistPainter.cxx.

◆ PaintErrors()

void THistPainter::PaintErrors ( Option_t option)
virtual

Draw 1D histograms error bars.

Definition at line 6228 of file THistPainter.cxx.

◆ PaintFrame()

void THistPainter::PaintFrame ( )
virtual

Calculate range and clear pad (canvas).

Definition at line 6682 of file THistPainter.cxx.

◆ PaintFunction()

void THistPainter::PaintFunction ( Option_t option)
virtual

[Paint functions associated to an histogram.](#HP28")

Definition at line 6705 of file THistPainter.cxx.

◆ PaintH3()

void THistPainter::PaintH3 ( Option_t option = "")
virtual

Control function to draw a 3D histograms.

Definition at line 6887 of file THistPainter.cxx.

◆ PaintH3Box()

void THistPainter::PaintH3Box ( Int_t  iopt)
virtual

Control function to draw a 3D histogram with boxes.

Definition at line 7376 of file THistPainter.cxx.

◆ PaintH3BoxRaster()

void THistPainter::PaintH3BoxRaster ( )
virtual

Control function to draw a 3D histogram with boxes.

Definition at line 7551 of file THistPainter.cxx.

◆ PaintH3Iso()

void THistPainter::PaintH3Iso ( )
virtual

Control function to draw a 3D histogram with Iso Surfaces.

Definition at line 7737 of file THistPainter.cxx.

◆ PaintHighlightBin()

void THistPainter::PaintHighlightBin ( Option_t option = "")
virtual

Paint highlight bin as TBox object.

Definition at line 3841 of file THistPainter.cxx.

◆ PaintHist()

void THistPainter::PaintHist ( Option_t option)
virtual

Control routine to draw 1D histograms

Definition at line 6745 of file THistPainter.cxx.

◆ PaintInit()

Int_t THistPainter::PaintInit ( )
virtual

Compute histogram parameters used by the drawing routines.

Definition at line 6975 of file THistPainter.cxx.

◆ PaintInitH()

Int_t THistPainter::PaintInitH ( )
virtual

Compute histogram parameters used by the drawing routines for a rotated pad.

Definition at line 7211 of file THistPainter.cxx.

◆ PaintLego()

void THistPainter::PaintLego ( Option_t option)
virtual

Control function to draw a 2D histogram as a lego plot.

Definition at line 7855 of file THistPainter.cxx.

◆ PaintLegoAxis()

void THistPainter::PaintLegoAxis ( TGaxis axis,
Double_t  ang 
)
virtual

Draw the axis for legos and surface plots.

Definition at line 8072 of file THistPainter.cxx.

◆ PaintPalette()

void THistPainter::PaintPalette ( )
virtual

Paint the color palette on the right side of the pad.

Definition at line 8256 of file THistPainter.cxx.

◆ PaintScatterPlot()

void THistPainter::PaintScatterPlot ( Option_t option)
virtual

Control function to draw a 2D histogram as a scatter plot.

Definition at line 8295 of file THistPainter.cxx.

◆ PaintSpecialObjects()

void THistPainter::PaintSpecialObjects ( const TObject obj,
Option_t option 
)
static

Static function to paint special objects like vectors and matrices.

This function is called via gROOT->ProcessLine to paint these objects without having a direct dependency of the graphics or histogramming system.

Definition at line 8411 of file THistPainter.cxx.

◆ PaintStat()

void THistPainter::PaintStat ( Int_t  dostat,
TF1 fit 
)
virtual

Draw the statistics box for 1D and profile histograms.

Implements TVirtualHistPainter.

Definition at line 8449 of file THistPainter.cxx.

◆ PaintStat2()

void THistPainter::PaintStat2 ( Int_t  dostat,
TF1 fit 
)
virtual

Draw the statistics box for 2D histograms.

Definition at line 8666 of file THistPainter.cxx.

◆ PaintStat3()

void THistPainter::PaintStat3 ( Int_t  dostat,
TF1 fit 
)
virtual

Draw the statistics box for 3D histograms.

Definition at line 8883 of file THistPainter.cxx.

◆ PaintSurface()

void THistPainter::PaintSurface ( Option_t option)
virtual

Control function to draw a 2D histogram as a surface plot.

Definition at line 9098 of file THistPainter.cxx.

◆ PaintTable()

void THistPainter::PaintTable ( Option_t option)
virtual

Control function to draw 2D/3D histograms (tables).

Definition at line 9475 of file THistPainter.cxx.

◆ PaintText()

void THistPainter::PaintText ( Option_t option)
virtual

Control function to draw a 1D/2D histograms with the bin values.

Definition at line 9906 of file THistPainter.cxx.

◆ PaintTF3()

void THistPainter::PaintTF3 ( )
virtual

Control function to draw a 3D implicit functions.

Definition at line 10002 of file THistPainter.cxx.

◆ PaintTH2PolyBins()

void THistPainter::PaintTH2PolyBins ( Option_t option)
virtual

Control function to draw a TH2Poly bins' contours.

  • option = "F" draw the bins as filled areas.
  • option = "L" draw the bins as line.
  • option = "P" draw the bins as markers.

Definition at line 9562 of file THistPainter.cxx.

◆ PaintTH2PolyColorLevels()

void THistPainter::PaintTH2PolyColorLevels ( Option_t option)
virtual

Control function to draw a TH2Poly as a color plot.

Definition at line 9638 of file THistPainter.cxx.

◆ PaintTH2PolyScatterPlot()

void THistPainter::PaintTH2PolyScatterPlot ( Option_t option)
virtual

Control function to draw a TH2Poly as a scatter plot.

Definition at line 9735 of file THistPainter.cxx.

◆ PaintTH2PolyText()

void THistPainter::PaintTH2PolyText ( Option_t option)
virtual

Control function to draw a TH2Poly as a text plot.

Definition at line 9848 of file THistPainter.cxx.

◆ PaintTitle()

void THistPainter::PaintTitle ( )
virtual

Draw the histogram title.

The title is drawn according to the title alignment returned by GetTitleAlign(). It is a 2 digits integer): hv

where h is the horizontal alignment and v is the vertical alignment.

  • h can get the values 1 2 3 for left, center, and right
  • v can get the values 1 2 3 for bottom, middle and top

for instance the default alignment is: 13 (left top)

Definition at line 10072 of file THistPainter.cxx.

◆ PaintTriangles()

void THistPainter::PaintTriangles ( Option_t option)
virtual

Control function to draw a table using Delaunay triangles.

Definition at line 9347 of file THistPainter.cxx.

◆ ProcessMessage()

void THistPainter::ProcessMessage ( const char *  mess,
const TObject obj 
)
virtual

Process message mess.

Implements TVirtualHistPainter.

Definition at line 10160 of file THistPainter.cxx.

◆ ProjectAitoff2xy()

Int_t THistPainter::ProjectAitoff2xy ( Double_t  l,
Double_t  b,
Double_t Al,
Double_t Ab 
)
static

Static function.

Convert Right Ascension, Declination to X,Y using an AITOFF projection. This procedure can be used to create an all-sky map in Galactic coordinates with an equal-area Aitoff projection. Output map coordinates are zero longitude centered. Also called Hammer-Aitoff projection (first presented by Ernst von Hammer in 1892)

source: GMT

code from Ernst-Jan Buis

Definition at line 10189 of file THistPainter.cxx.

◆ ProjectMercator2xy()

Int_t THistPainter::ProjectMercator2xy ( Double_t  l,
Double_t  b,
Double_t Al,
Double_t Ab 
)
static

Static function.

Probably the most famous of the various map projections, the Mercator projection takes its name from Mercator who presented it in 1569. It is a cylindrical, conformal projection with no distortion along the equator. The Mercator projection has been used extensively for world maps in which the distortion towards the polar regions grows rather large, thus incorrectly giving the impression that, for example, Greenland is larger than South America. In reality, the latter is about eight times the size of Greenland. Also, the Former Soviet Union looks much bigger than Africa or South America. One may wonder whether this illusion has had any influence on U.S. foreign policy.' (Source: GMT) code from Ernst-Jan Buis

Definition at line 10224 of file THistPainter.cxx.

◆ ProjectParabolic2xy()

Int_t THistPainter::ProjectParabolic2xy ( Double_t  l,
Double_t  b,
Double_t Al,
Double_t Ab 
)
static

Static function code from Ernst-Jan Buis.

Definition at line 10247 of file THistPainter.cxx.

◆ ProjectSinusoidal2xy()

Int_t THistPainter::ProjectSinusoidal2xy ( Double_t  l,
Double_t  b,
Double_t Al,
Double_t Ab 
)
static

Static function code from Ernst-Jan Buis.

Definition at line 10236 of file THistPainter.cxx.

◆ RecalculateRange()

void THistPainter::RecalculateRange ( )
virtual

Recompute the histogram range following graphics operations.

Definition at line 10258 of file THistPainter.cxx.

◆ RecursiveRemove()

virtual void THistPainter::RecursiveRemove ( TObject obj)
inlinevirtual

Recursively remove this object from a list.

Typically implemented by classes that can contain multiple references to a same object.

Reimplemented from TObject.

Definition at line 139 of file THistPainter.h.

◆ SetHighlight()

void THistPainter::SetHighlight ( )
virtual

Set highlight (enable/disable) mode for fH.

Implements TVirtualHistPainter.

Definition at line 3790 of file THistPainter.cxx.

◆ SetHistogram()

void THistPainter::SetHistogram ( TH1 h)
virtual

Set current histogram to h

Implements TVirtualHistPainter.

Definition at line 10369 of file THistPainter.cxx.

◆ SetShowProjection()

void THistPainter::SetShowProjection ( const char *  option,
Int_t  nbins 
)
virtual

Set projection.

Implements TVirtualHistPainter.

Definition at line 10613 of file THistPainter.cxx.

◆ SetStack()

virtual void THistPainter::SetStack ( TList stack)
inlinevirtual

Implements TVirtualHistPainter.

Definition at line 142 of file THistPainter.h.

◆ ShowProjection3()

void THistPainter::ShowProjection3 ( Int_t  px,
Int_t  py 
)
virtual

Show projection (specified by fShowProjection) of a TH3.

The drawing option for the projection is in fShowOption.

First implementation; R.Brun

Full implementation: Tim Tran (timtr.nosp@m.an@j.nosp@m.lab.o.nosp@m.rg) April 2006

Definition at line 10811 of file THistPainter.cxx.

◆ ShowProjectionX()

void THistPainter::ShowProjectionX ( Int_t  px,
Int_t  py 
)
virtual

Show projection onto X.

Definition at line 10640 of file THistPainter.cxx.

◆ ShowProjectionY()

void THistPainter::ShowProjectionY ( Int_t  px,
Int_t  py 
)
virtual

Show projection onto Y.

Definition at line 10723 of file THistPainter.cxx.

◆ TableInit()

Int_t THistPainter::TableInit ( )
virtual

Initialize various options to draw 2D histograms.

Definition at line 10383 of file THistPainter.cxx.

Member Data Documentation

◆ fCuts

TCutG* THistPainter::fCuts[kMaxCuts]
protected

Definition at line 62 of file THistPainter.h.

◆ fCutsOpt

Int_t THistPainter::fCutsOpt[kMaxCuts]
protected

Definition at line 61 of file THistPainter.h.

◆ fFunctions

TList* THistPainter::fFunctions
protected

Definition at line 54 of file THistPainter.h.

◆ fGraph2DPainter

TGraph2DPainter* THistPainter::fGraph2DPainter
protected

Definition at line 56 of file THistPainter.h.

◆ fH

TH1* THistPainter::fH
protected

Definition at line 50 of file THistPainter.h.

◆ fLego

TPainter3dAlgorithms* THistPainter::fLego
protected

Definition at line 55 of file THistPainter.h.

◆ fNcuts

Int_t THistPainter::fNcuts
protected

Definition at line 60 of file THistPainter.h.

◆ fObjectInfo

TString THistPainter::fObjectInfo
mutableprivate

Definition at line 70 of file THistPainter.h.

◆ fPie

TPie* THistPainter::fPie
protected

Definition at line 57 of file THistPainter.h.

◆ fShowOption

TString THistPainter::fShowOption
protected

Definition at line 65 of file THistPainter.h.

◆ fShowProjection

Int_t THistPainter::fShowProjection
protected

Definition at line 64 of file THistPainter.h.

◆ fStack

TList* THistPainter::fStack
protected

Definition at line 63 of file THistPainter.h.

◆ fXaxis

TAxis* THistPainter::fXaxis
protected

Definition at line 51 of file THistPainter.h.

◆ fXbuf

Double_t* THistPainter::fXbuf
protected

Definition at line 58 of file THistPainter.h.

◆ fXHighlightBin

Int_t THistPainter::fXHighlightBin
protected

Definition at line 66 of file THistPainter.h.

◆ fYaxis

TAxis* THistPainter::fYaxis
protected

Definition at line 52 of file THistPainter.h.

◆ fYbuf

Double_t* THistPainter::fYbuf
protected

Definition at line 59 of file THistPainter.h.

◆ fYHighlightBin

Int_t THistPainter::fYHighlightBin
protected

Definition at line 67 of file THistPainter.h.

◆ fZaxis

TAxis* THistPainter::fZaxis
protected

Definition at line 53 of file THistPainter.h.

Libraries for THistPainter:
[legend]

The documentation for this class was generated from the following files: