Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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
12void candleplotwhiskers() {
13 auto c1 = new TCanvas("c1","Candle Presets",700,800);
14 c1->Divide(1,2);
15
16 auto rng = new TRandom();
17 auto h1 = new TH2I("h1","Gaus",100,-5,5,1,0,1);
18 auto 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 = rng->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 auto mygaus_1_middle = new TF1("mygaus_1_middle","gaus",q[0],q[2]);
40 auto mygaus_1_left = new TF1("mygaus_1_left","gaus",q[0]-1.5*iqr,q[0]);
41 mygaus_1_left->SetLineColor(kGreen);
42 auto 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->SetFillStyle(0);
47 h1->Draw("candley2 scat");
48
49 c1->cd(2);
50 h2->Draw("");
51 h2->Fit("mygaus_1_left","R");
52 mygaus_1_left->Draw("same");
53 auto l3 = new TLine(q[0]-1.5*iqr,0,q[0]-1.5*iqr,mygaus_1_left->Eval(q[0]-1.5*iqr));
54 l3->SetLineColor(kGreen); l3->SetLineWidth(2); l3->Draw("");
55 auto l1 = new TLine(q[0] ,0,q[0] ,mygaus_1_left->Eval(q[0]));
56 l1->SetLineWidth(2); l1->SetLineColor(kGreen); l1->Draw("");
57
58 h2->Fit("mygaus_1_right","R","");
59 mygaus_1_right->Draw("same");
60 auto l4 = new TLine(q[2]+1.5*iqr,0,q[2]+1.5*iqr,mygaus_1_left->Eval(q[2]+1.5*iqr));
61 l4->SetLineColor(kGreen); l4->SetLineWidth(2); l4->Draw("");
62 auto l5 = new TLine(q[2] ,0,q[2] ,mygaus_1_right->Eval(q[2]));
63 l5->SetLineWidth(2); l5->SetLineColor(kGreen); l5->Draw("");
64
65 h2->Fit("mygaus_1_middle","R");
66 mygaus_1_middle->Draw("same");
67
68 //In principal one could calculate these values by h2->Integral() as well
69 TText t;
70 t.SetTextFont(42);
71 t.DrawText(0,mygaus_1_middle->Eval(0)/2,"50%");
72 t.DrawText(-1.5,mygaus_1_middle->Eval(-1.5)/2,"24.65%");
73 t.DrawText(+1,mygaus_1_middle->Eval(+1.5)/2,"24.65%");
74 t.DrawText(q[0]-1.5*iqr,1000,Form("%.3f",q[0]-1.5*iqr))->SetTextAngle(90);
75 t.DrawText(q[2]+1.5*iqr,1000,Form("%.3f",q[2]+1.5*iqr))->SetTextAngle(90);
76 t.DrawText(q[0],1000,Form("%.3f",q[0]))->SetTextAngle(90);
77 t.DrawText(q[2],1000,Form("%.3f",q[2]))->SetTextAngle(90);
78}
double Double_t
Definition RtypesCore.h:59
@ kGreen
Definition Rtypes.h:66
winID h TVirtualViewer3D TVirtualGLPainter p
float * q
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2489
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
virtual void SetTextAngle(Float_t tangle=0)
Set the text angle.
Definition TAttText.h:43
virtual void SetTextFont(Font_t tfont=62)
Set the text font.
Definition TAttText.h:46
The Canvas class.
Definition TCanvas.h:23
1-Dim function class
Definition TF1.h:233
1-D histogram with an int per channel (see TH1 documentation)
Definition TH1.h:541
TAxis * GetXaxis()
Definition TH1.h:325
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition TH1.cxx:3346
void Draw(Option_t *option="") override
Draw this histogram with options.
Definition TH1.cxx:3068
2-D histogram with an int per channel (see TH1 documentation)
Definition TH2.h:226
Use the TLine constructor to create a simple line.
Definition TLine.h:22
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition TNamed.cxx:164
This is the base class for the ROOT Random number generators.
Definition TRandom.h:27
Base class for several text objects.
Definition TText.h:22
virtual TText * DrawText(Double_t x, Double_t y, const char *text)
Draw this text with new coordinates.
Definition TText.cxx:176
return c1
Definition legend1.C:41
TH1F * h1
Definition legend1.C:5