Logo ROOT   6.08/07
Reference Guide
candleplotwhiskers.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_hist
3 /// \notebook
4 /// Example of candle plot showing the whiskers definition.
5 ///
6 /// \macro_image
7 /// \macro_output
8 /// \macro_code
9 ///
10 /// \author Georg Troska
11 
12 void candleplotwhiskers() {
13  TCanvas *c1 = new TCanvas("c1","Candle Presets",700,800);
14  c1->Divide(1,2);
15 
16  TRandom *randnum = new TRandom();
17  TH2I *h1 = new TH2I("h1","Gaus",100,-5,5,1,0,1);
18  TH1I *h2 = new TH1I("h2","Gaus",100,-5,5);
19 
20  h1->GetXaxis()->SetTitle("Standard deviation #sigma");
21  h2->GetXaxis()->SetTitle("Standard deviation #sigma");
22  h2->GetYaxis()->SetTitle("dN/d#sigma");
23 
24  float myRand;
25  for (int i = 0; i < 100000; i++) {
26  myRand = randnum->Gaus(0,1);
27  h1->Fill(myRand,0);
28  h2->Fill(myRand);
29  }
30 
31  Double_t *q = new Double_t[3];
32  Double_t *p = new Double_t[3];
33  q[0] = 0.; q[1] = 0.; q[2] = 0.;
34  p[0] = 0.25; p[1] = 0.5; p[2] = 0.75;
35 
36  h2->GetQuantiles(3,q,p);
37  cout << "Q1 (-25%): " << q[0] << " Median: " << q[1] << " Q3 (+25%): " << q[2] << endl;
38  double iqr = q[2]-q[0];
39  TF1 *mygaus_1_middle = new TF1("mygaus_1_middle","gaus",q[0],q[2]);
40  TF1 *mygaus_1_left = new TF1("mygaus_1_left","gaus",q[0]-1.5*iqr,q[0]);
41  mygaus_1_left->SetLineColor(kGreen);
42  TF1 *mygaus_1_right = new TF1("mygaus_1_right","gaus",q[2],q[2]+1.5*iqr);
43  mygaus_1_right->SetLineColor(kGreen);
44  c1->cd(1);
45  h1->SetLineWidth(3);
46  h1->Draw("candley2 scat");
47 
48  c1->cd(2);
49  h2->Draw("");
50  h2->Fit("mygaus_1_left","R");
51  mygaus_1_left->Draw("same");
52  TLine *l3 = new TLine(q[0]-1.5*iqr,0,q[0]-1.5*iqr,mygaus_1_left->Eval(q[0]-1.5*iqr));
53  l3->SetLineColor(kGreen); l3->SetLineWidth(2); l3->Draw("");
54  TLine *l1 = new TLine(q[0] ,0,q[0] ,mygaus_1_left->Eval(q[0]));
55  l1->SetLineWidth(2); l1->SetLineColor(kGreen); l1->Draw("");
56 
57  h2->Fit("mygaus_1_right","R","");
58  mygaus_1_right->Draw("same");
59  TLine *l4 = new TLine(q[2]+1.5*iqr,0,q[2]+1.5*iqr,mygaus_1_left->Eval(q[2]+1.5*iqr));
60  l4->SetLineColor(kGreen); l4->SetLineWidth(2); l4->Draw("");
61  TLine *l5 = new TLine(q[2] ,0,q[2] ,mygaus_1_right->Eval(q[2]));
62  l5->SetLineWidth(2); l5->SetLineColor(kGreen); l5->Draw("");
63 
64  h2->Fit("mygaus_1_middle","R");
65  mygaus_1_middle->Draw("same");
66 
67  //In principal one could calculate these values by h2->Integral() as well
68  TText *t = new TText(); t->SetTextFont(42);
69  t->DrawText(0,mygaus_1_middle->Eval(0)/2,"50%");
70  t->DrawText(-1.5,mygaus_1_middle->Eval(-1.5)/2,"24.65%");
71  t->DrawText(+1,mygaus_1_middle->Eval(+1.5)/2,"24.65%");
72  t->DrawText(q[0]-1.5*iqr,1000,Form("%.3f",q[0]-1.5*iqr))->SetTextAngle(90);
73  t->DrawText(q[2]+1.5*iqr,1000,Form("%.3f",q[2]+1.5*iqr))->SetTextAngle(90);
74  t->DrawText(q[0],1000,Form("%.3f",q[0]))->SetTextAngle(90);
75  t->DrawText(q[2],1000,Form("%.3f",q[2]))->SetTextAngle(90);
76 }
virtual void SetLineWidth(Width_t lwidth)
Set the line width.
Definition: TAttLine.h:49
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition: TH1.cxx:3125
THist< 1, int, THistStatContent > TH1I
Definition: THist.hxx:304
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:235
return c1
Definition: legend1.C:41
TVirtualPad * cd(Int_t subpadnumber=0)
Set current canvas & pad.
Definition: TCanvas.cxx:659
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:4190
virtual void Draw(Option_t *option="")
Draw this function with its current attributes.
Definition: TF1.cxx:1086
virtual void Draw(Option_t *option="")
Default Draw method for all objects.
Definition: TObject.cxx:255
Definition: Rtypes.h:61
virtual void SetTextFont(Font_t tfont=62)
Set the text font.
Definition: TAttText.h:51
This is the base class for the ROOT Random number generators.
Definition: TRandom.h:31
Base class for several text objects.
Definition: TText.h:33
TH1F * h1
Definition: legend1.C:5
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition: TAttLine.h:46
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2851
tomato 1-D histogram with an int per channel (see TH1 documentation)}
Definition: TH1.h:534
virtual void SetTextAngle(Float_t tangle=0)
Set the text angle.
Definition: TAttText.h:48
char * Form(const char *fmt,...)
A simple line.
Definition: TLine.h:33
TAxis * GetYaxis()
Definition: TH1.h:325
The Canvas class.
Definition: TCanvas.h:41
virtual Double_t Eval(Double_t x, Double_t y=0, Double_t z=0, Double_t t=0) const
Evaluate this function.
Definition: TF1.cxx:1196
double Double_t
Definition: RtypesCore.h:55
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:1089
1-Dim function class
Definition: TF1.h:149
virtual TText * DrawText(Double_t x, Double_t y, const char *text)
Draw this text with new coordinates.
Definition: TText.cxx:175
THist< 2, int, THistStatContent > TH2I
Definition: THist.hxx:310
tomato 2-D histogram with an int per channel (see TH1 documentation)}
Definition: TH2.h:216
float * q
Definition: THbookFile.cxx:87
Int_t Fill(Double_t)
Invalid Fill method.
Definition: TH2.cxx:292
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:3563
TAxis * GetXaxis()
Definition: TH1.h:324