Logo ROOT  
Reference Guide
quantiles.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_math
3/// \notebook -js
4/// Demo for quantiles
5///
6/// \macro_image
7/// \macro_code
8///
9/// \authors Rene Brun, Eddy Offermann
10
11void quantiles() {
12 const Int_t nq = 100;
13 const Int_t nshots = 10;
14 Double_t xq[nq]; // position where to compute the quantiles in [0,1]
15 Double_t yq[nq]; // array to contain the quantiles
16 for (Int_t i=0;i<nq;i++) xq[i] = Float_t(i+1)/nq;
17
18 TGraph *gr70 = new TGraph(nshots);
19 TGraph *gr90 = new TGraph(nshots);
20 TGraph *gr98 = new TGraph(nshots);
21 TH1F *h = new TH1F("h","demo quantiles",50,-3,3);
22
23 for (Int_t shot=0;shot<nshots;shot++) {
24 h->FillRandom("gaus",50);
25 h->GetQuantiles(nq,yq,xq);
26 gr70->SetPoint(shot,shot+1,yq[70]);
27 gr90->SetPoint(shot,shot+1,yq[90]);
28 gr98->SetPoint(shot,shot+1,yq[98]);
29 }
30
31 //show the original histogram in the top pad
32 TCanvas *c1 = new TCanvas("c1","demo quantiles",10,10,600,900);
33 c1->Divide(1,3);
34 c1->cd(1);
35 h->SetFillColor(38);
36 h->Draw();
37
38 // show the final quantiles in the middle pad
39 c1->cd(2);
40 gPad->SetGrid();
41 TGraph *gr = new TGraph(nq,xq,yq);
42 gr->SetTitle("final quantiles");
43 gr->SetMarkerStyle(21);
45 gr->SetMarkerSize(0.3);
46 gr->Draw("ap");
47
48 // show the evolution of some quantiles in the bottom pad
49 c1->cd(3);
50 gPad->DrawFrame(0,0,nshots+1,3.2);
51 gPad->SetGrid();
52 gr98->SetMarkerStyle(22);
53 gr98->SetMarkerColor(kRed);
54 gr98->Draw("lp");
55 gr90->SetMarkerStyle(21);
56 gr90->SetMarkerColor(kBlue);
57 gr90->Draw("lp");
58 gr70->SetMarkerStyle(20);
60 gr70->Draw("lp");
61 // add a legend
62 TLegend *legend = new TLegend(0.85,0.74,0.95,0.95);
63 legend->SetTextFont(72);
64 legend->SetTextSize(0.05);
65 legend->AddEntry(gr98," q98","lp");
66 legend->AddEntry(gr90," q90","lp");
67 legend->AddEntry(gr70," q70","lp");
68 legend->Draw();
69}
#define h(i)
Definition: RSha256.hxx:106
int Int_t
Definition: RtypesCore.h:43
double Double_t
Definition: RtypesCore.h:57
float Float_t
Definition: RtypesCore.h:55
@ kRed
Definition: Rtypes.h:64
@ kMagenta
Definition: Rtypes.h:64
@ kBlue
Definition: Rtypes.h:64
#define gPad
Definition: TVirtualPad.h:287
virtual void SetMarkerColor(Color_t mcolor=1)
Set the marker color.
Definition: TAttMarker.h:38
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
Definition: TAttMarker.h:40
virtual void SetMarkerSize(Size_t msize=1)
Set the marker size.
Definition: TAttMarker.h:41
virtual void SetTextFont(Font_t tfont=62)
Set the text font.
Definition: TAttText.h:45
virtual void SetTextSize(Float_t tsize=1)
Set the text size.
Definition: TAttText.h:46
The Canvas class.
Definition: TCanvas.h:27
A TGraph is an object made of two arrays X and Y with npoints each.
Definition: TGraph.h:41
virtual void SetPoint(Int_t i, Double_t x, Double_t y)
Set x and y values for point number i.
Definition: TGraph.cxx:2269
virtual void SetTitle(const char *title="")
Change (i.e.
Definition: TGraph.cxx:2324
virtual void Draw(Option_t *chopt="")
Draw this graph with its current attributes.
Definition: TGraph.cxx:760
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:571
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
return c1
Definition: legend1.C:41
TGraphErrors * gr
Definition: legend1.C:25