Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
hlquantiles.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_math
3/// Demo for quantiles (with highlight mode)
4///
5/// \macro_image
6/// \macro_code
7///
8/// \authors Rene Brun, Eddy Offermann, Jan Musinsky
9
10TList *lq = nullptr;
11TGraph *gr = nullptr;
12
13void HighlightQuantile(TVirtualPad *pad, TObject *obj, Int_t ihp, Int_t y)
14{
15 // show the evolution of all quantiles in the bottom pad
16 if (obj != gr) return;
17 if (ihp == -1) return;
18
19 TVirtualPad *savepad = gPad;
20 pad->GetCanvas()->cd(3);
21 lq->At(ihp)->Draw("alp");
22 gPad->Update();
23 if (savepad) savepad->cd();
24}
25
26
27void hlquantiles() {
28 const Int_t nq = 100;
29 const Int_t nshots = 10;
30 Double_t xq[nq]; // position where to compute the quantiles in [0,1]
31 Double_t yq[nq]; // array to contain the quantiles
32 for (Int_t i=0;i<nq;i++) xq[i] = Float_t(i+1)/nq;
33
34 TGraph *gr70 = new TGraph(nshots);
35 TGraph *gr90 = new TGraph(nshots);
36 TGraph *gr98 = new TGraph(nshots);
37 TGraph *grq[nq];
38 for (Int_t ig = 0; ig < nq; ig++)
39 grq[ig] = new TGraph(nshots);
40 TH1F *h = new TH1F("h","demo quantiles",50,-3,3);
41
42 for (Int_t shot=0;shot<nshots;shot++) {
43 h->FillRandom("gaus",50);
44 h->GetQuantiles(nq,yq,xq);
45 gr70->SetPoint(shot,shot+1,yq[70]);
46 gr90->SetPoint(shot,shot+1,yq[90]);
47 gr98->SetPoint(shot,shot+1,yq[98]);
48 for (Int_t ig = 0; ig < nq; ig++)
49 grq[ig]->SetPoint(shot,shot+1,yq[ig]);
50 }
51
52 //show the original histogram in the top pad
53 TCanvas *c1 = new TCanvas("c1","demo quantiles",10,10,600,900);
54 c1->HighlightConnect("HighlightQuantile(TVirtualPad*,TObject*,Int_t,Int_t)");
55 c1->SetFillColor(41);
56 c1->Divide(1,3);
57 c1->cd(1);
58 h->SetFillColor(38);
59 h->Draw();
60
61 // show the final quantiles in the middle pad
62 c1->cd(2);
63 gPad->SetFrameFillColor(33);
64 gPad->SetGrid();
65 gr = new TGraph(nq,xq,yq);
66 gr->SetTitle("final quantiles");
67 gr->SetMarkerStyle(21);
69 gr->SetMarkerSize(0.3);
70 gr->Draw("ap");
71
72 // prepare quantiles
73 lq = new TList();
74 for (Int_t ig = 0; ig < nq; ig++) {
75 grq[ig]->SetMinimum(gr->GetYaxis()->GetXmin());
76 grq[ig]->SetMaximum(gr->GetYaxis()->GetXmax());
77 grq[ig]->SetMarkerStyle(23);
78 grq[ig]->SetMarkerColor(ig%100);
79 grq[ig]->SetTitle(TString::Format("q%02d", ig));
80 lq->Add(grq[ig]);
81 }
82
83 TText *info = new TText(0.1, 2.4, "please move the mouse over the graph");
84 info->SetTextSize(0.08);
86 info->SetBit(kCannotPick);
87 info->Draw();
88
90
91 // show the evolution of some quantiles in the bottom pad
92 c1->cd(3);
93 gPad->SetFrameFillColor(17);
94 gPad->DrawFrame(0,0,nshots+1,3.2);
95 gPad->SetGrid();
96 gr98->SetMarkerStyle(22);
97 gr98->SetMarkerColor(kRed);
98 gr98->Draw("lp");
99 gr90->SetMarkerStyle(21);
100 gr90->SetMarkerColor(kBlue);
101 gr90->Draw("lp");
102 gr70->SetMarkerStyle(20);
104 gr70->Draw("lp");
105 // add a legend
106 TLegend *legend = new TLegend(0.85,0.74,0.95,0.95);
107 legend->SetTextFont(72);
108 legend->SetTextSize(0.05);
109 legend->AddEntry(gr98," q98","lp");
110 legend->AddEntry(gr90," q90","lp");
111 legend->AddEntry(gr70," q70","lp");
112 legend->Draw();
113}
114
#define h(i)
Definition RSha256.hxx:106
int Int_t
Definition RtypesCore.h:45
float Float_t
Definition RtypesCore.h:57
double Double_t
Definition RtypesCore.h:59
@ kRed
Definition Rtypes.h:66
@ kMagenta
Definition Rtypes.h:66
@ kBlue
Definition Rtypes.h:66
int * lq
@ kCannotPick
Definition TObject.h:374
#define gPad
virtual void SetMarkerColor(Color_t mcolor=1)
Set the marker color.
Definition TAttMarker.h:38
virtual Color_t GetMarkerColor() const
Return the marker color.
Definition TAttMarker.h:31
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:45
virtual void SetTextColor(Color_t tcolor=1)
Set the text color.
Definition TAttText.h:44
virtual void SetTextFont(Font_t tfont=62)
Set the text font.
Definition TAttText.h:46
virtual void SetTextSize(Float_t tsize=1)
Set the text size.
Definition TAttText.h:47
Double_t GetXmax() const
Definition TAxis.h:140
Double_t GetXmin() const
Definition TAxis.h:139
The Canvas class.
Definition TCanvas.h:23
TVirtualPad * cd(Int_t subpadnumber=0) override
Set current canvas & pad.
Definition TCanvas.cxx:716
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:2319
virtual void SetMaximum(Double_t maximum=-1111)
Set the maximum of the graph.
Definition TGraph.cxx:2301
virtual void SetHighlight(Bool_t set=kTRUE)
Set highlight (enable/disable) mode for the graph by default highlight mode is disable.
Definition TGraph.cxx:2288
void Draw(Option_t *chopt="") override
Draw this graph with its current attributes.
Definition TGraph.cxx:809
TAxis * GetYaxis() const
Get y axis of the graph.
Definition TGraph.cxx:1553
void SetTitle(const char *title="") override
Change (i.e.
Definition TGraph.cxx:2374
virtual void SetMinimum(Double_t minimum=-1111)
Set the minimum of the graph.
Definition TGraph.cxx:2310
1-D histogram with a float per channel (see TH1 documentation)
Definition TH1.h:621
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:317
void Draw(Option_t *option="") override
Draw this legend with its current attributes.
Definition TLegend.cxx:422
A doubly linked list.
Definition TList.h:38
Mother of all ROOT objects.
Definition TObject.h:41
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:780
virtual void Draw(Option_t *option="")
Default Draw method for all objects.
Definition TObject.cxx:274
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:2378
Base class for several text objects.
Definition TText.h:22
TVirtualPad is an abstract base class for the Pad and Canvas classes.
Definition TVirtualPad.h:51
virtual TVirtualPad * cd(Int_t subpadnumber=0)=0
virtual TCanvas * GetCanvas() const =0
Double_t y[n]
Definition legend1.C:17
return c1
Definition legend1.C:41
TGraphErrors * gr
Definition legend1.C:25