Logo ROOT   6.10/09
Reference Guide
ratioplot2.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_hist
3 /// \notebook
4 /// Example of a fit residual plot.
5 ///
6 /// Creates a histogram filled with random numbers from a gaussian distribution
7 /// and fits it with a standard gaussian function. The result is passed to the `TRatioPlot`
8 /// constructor. Additionally, after calling `TRatioPlot::Draw` the upper and lower y axis
9 /// titles are modified.
10 /// Confidence interval bands are automatically drawn on the bottom (but can be disabled by draw option `nobands`.
11 ///
12 /// \macro_image
13 /// \macro_code
14 ///
15 /// \author Paul Gessinger
16 
17 void ratioplot2() {
18  gStyle->SetOptStat(0);
19  auto c1 = new TCanvas("c1", "fit residual simple");
20  auto h1 = new TH1D("h1", "h1", 50, -5, 5);
21  h1->FillRandom("gaus", 2000);
22  h1->Fit("gaus");
23  h1->GetXaxis()->SetTitle("x");
24  c1->Clear(); // Fit does not draw into correct pad
25  auto rp1 = new TRatioPlot(h1);
26  rp1->Draw();
27  rp1->GetLowerRefYaxis()->SetTitle("ratio");
28  rp1->GetUpperRefYaxis()->SetTitle("entries");
29  c1->Update();
30 }
return c1
Definition: legend1.C:41
R__EXTERN TStyle * gStyle
Definition: TStyle.h:402
TH1F * h1
Definition: legend1.C:5
virtual void FillRandom(const char *fname, Int_t ntimes=5000)
Fill histogram following distribution in function fname.
Definition: TH1.cxx:3294
The Canvas class.
Definition: TCanvas.h:31
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:1267
THist< 1, double, THistStatContent, THistStatUncertainty > TH1D
Definition: THist.hxx:310
Class for displaying ratios, differences and fit residuals.
Definition: TRatioPlot.h:44
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition: TNamed.cxx:155
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:3564
TAxis * GetXaxis()
Definition: TH1.h:300