ROOT
Version v6.34
master
v6.32
Reference Guide
▼
ROOT
ROOT Reference Documentation
Tutorials
▼
Functional Parts
►
Core ROOT classes
►
std Extension classes
►
Parallelized classes
►
The Geometry Package
►
Graphics
►
Event display with ROOT7
►
GUI
►
Web Widgets
►
Web Display
►
Histogram Library
►
Input/Output Library
►
Math
►
N-D parametric functions
►
VecOps
►
Monte Carlo
►
HTTP server
►
PROOF
►
TMVA
►
RooFit
►
Dataframe
►
ROOT7 classes
►
NTuple-related classes
►
Tree Library
►
TreePlayer Library
▼
Tutorials
►
Histograms tutorials
►
Tree tutorials
►
Dataframe tutorials
►
ROOT 7 tutorials
►
FOAM tutorials
►
Containers tutorials
►
Event display tutorials
►
Event display ROOT7 tutorials
►
Geometry tutorials
►
Fast Fourier Transforms tutorials
►
Fit Tutorials
►
RooFit Tutorials
►
Graphs tutorials
►
Graphics tutorials
►
OpenGL tutorials
►
Tutorials specific to Mac/Cocoa
►
GUI tutorials
►
HistFactory Tutorials
►
HTTP tutorials
►
Image tutorials
►
IO tutorials
▼
Math tutorials
Bessel.C
Bessel.py
binomial.C
BreitWigner.C
ChebyshevPol.C
chi2test.C
CrystalBall.C
exampleFunction.py
exampleFunctor.C
exampleMultiRoot.C
exampleTKDE.C
FeldmanCousins.C
GammaFun.C
goftest.C
hlquantiles.C
kdTreeBinning.C
Legendre.C
Legendre.py
LegendreAssoc.C
limit.C
mathBeta.C
mathcoreCDF.C
mathcoreGenVector.C
mathcoreSpecFunc.C
mathcoreStatFunc.C
mathcoreStatFunc.py
mathcoreVectorCollection.C
mathcoreVectorFloatIO.C
mathcoreVectorIO.C
mathGammaNormal.C
mathLaplace.C
mathmoreIntegration.C
mathmoreIntegrationMultidim.C
mathStudent.C
multidimSampling.C
multivarGaus.C
normalDist.C
normalDist.py
permute.C
principal.C
principal.py
quantiles.C
quasirandom.C
Rolke.C
testrandom.C
tStudent.C
tStudent.py
TSVDUnfoldExample.C
vavilov.C
►
Matrix tutorials
►
Monte Carlo tutorials
►
Multicore tutorials
►
Net tutorials
►
Physics tutorials
►
PyRoot tutorials
►
Pythia tutorials
►
Quadratic programming package.
►
R tutorials
►
RooStats Tutorials
►
Spectrum tutorials
►
TSPlot tutorials
►
SQL tutorials
►
TMVA tutorials
►
TUnfold tutorials
►
Unuran tutorials
►
VecOps tutorials
►
FITS files interface tutorials
►
XML tutorials
►
Proof tutorials
►
TWebCanvas tutorials
►
Webgui tutorials
►
Legacy tutorials
demos.C
demoshelp.C
hsimple.C
rootlogoff.C
rootlogon.C
►
R Interface for Statistical Computing
►
Namespaces
►
All Classes
►
Files
Release Notes
•
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Properties
Friends
Macros
Modules
Pages
Loading...
Searching...
No Matches
hlquantiles.C File Reference
Tutorials
»
Math tutorials
Detailed Description
Demo for quantiles (with highlight mode)
TList
*
lq
=
nullptr
;
TGraph
*
gr
=
nullptr
;
void
HighlightQuantile
(
TVirtualPad
*
pad
,
TObject
*obj,
Int_t
ihp
,
Int_t
y
)
{
// show the evolution of all quantiles in the bottom pad
if
(obj !=
gr
)
return
;
if
(
ihp
== -1)
return
;
TVirtualPad
*
savepad
=
gPad
;
pad
->GetCanvas()->cd(3);
lq
->At(
ihp
)->Draw(
"alp"
);
gPad
->Update();
if
(
savepad
)
savepad
->cd();
}
void
hlquantiles
() {
const
Int_t
nq
= 100;
const
Int_t
nshots
= 10;
Double_t
xq
[
nq
];
// position where to compute the quantiles in [0,1]
Double_t
yq
[
nq
];
// array to contain the quantiles
for
(
Int_t
i=0;i<
nq
;i++)
xq
[i] =
Float_t
(i+1)/
nq
;
TGraph
*
gr70
=
new
TGraph
(
nshots
);
TGraph
*
gr90
=
new
TGraph
(
nshots
);
TGraph
*
gr98
=
new
TGraph
(
nshots
);
TGraph
*
grq
[
nq
];
for
(
Int_t
ig
= 0;
ig
<
nq
;
ig
++)
grq
[
ig
] =
new
TGraph
(
nshots
);
TH1F
*
h
=
new
TH1F
(
"h"
,
"demo quantiles"
,50,-3,3);
for
(
Int_t
shot
=0;
shot
<
nshots
;
shot
++) {
h
->FillRandom(
"gaus"
,50);
h
->GetQuantiles(
nq
,
yq
,
xq
);
gr70
->SetPoint(
shot
,
shot
+1,
yq
[70]);
gr90
->SetPoint(
shot
,
shot
+1,
yq
[90]);
gr98
->SetPoint(
shot
,
shot
+1,
yq
[98]);
for
(
Int_t
ig
= 0;
ig
<
nq
;
ig
++)
grq
[
ig
]->SetPoint(
shot
,
shot
+1,
yq
[
ig
]);
}
//show the original histogram in the top pad
TCanvas
*
c1
=
new
TCanvas
(
"c1"
,
"demo quantiles"
,10,10,600,900);
c1
->HighlightConnect(
"HighlightQuantile(TVirtualPad*,TObject*,Int_t,Int_t)"
);
c1
->SetFillColor(41);
c1
->Divide(1,3);
c1
->cd(1);
h
->SetFillColor(38);
h
->Draw();
// show the final quantiles in the middle pad
c1
->cd(2);
gPad
->SetFrameFillColor(33);
gPad
->SetGrid();
gr
=
new
TGraph
(
nq
,
xq
,
yq
);
gr
->
SetTitle
(
"final quantiles"
);
gr
->
SetMarkerStyle
(21);
gr
->
SetMarkerColor
(
kRed
);
gr
->
SetMarkerSize
(0.3);
gr
->
Draw
(
"ap"
);
// prepare quantiles
lq
=
new
TList
();
for
(
Int_t
ig
= 0;
ig
<
nq
;
ig
++) {
grq
[
ig
]->SetMinimum(
gr
->
GetYaxis
()->
GetXmin
());
grq
[
ig
]->SetMaximum(
gr
->
GetYaxis
()->
GetXmax
());
grq
[
ig
]->SetMarkerStyle(23);
grq
[
ig
]->SetMarkerColor(
ig
%100);
grq
[
ig
]->SetTitle(
TString::Format
(
"q%02d"
,
ig
));
lq
->Add(
grq
[
ig
]);
}
TText
*
info
=
new
TText
(0.1, 2.4,
"please move the mouse over the graph"
);
info
->SetTextSize(0.08);
info
->SetTextColor(
gr
->
GetMarkerColor
());
info
->SetBit(
kCannotPick
);
info
->Draw();
gr
->
SetHighlight
();
// show the evolution of some quantiles in the bottom pad
c1
->cd(3);
gPad
->SetFrameFillColor(17);
gPad
->DrawFrame(0,0,
nshots
+1,3.2);
gPad
->SetGrid();
gr98
->SetMarkerStyle(22);
gr98
->SetMarkerColor(
kRed
);
gr98
->Draw(
"lp"
);
gr90
->SetMarkerStyle(21);
gr90
->SetMarkerColor(
kBlue
);
gr90
->Draw(
"lp"
);
gr70
->SetMarkerStyle(20);
gr70
->SetMarkerColor(
kMagenta
);
gr70
->Draw(
"lp"
);
// add a legend
TLegend
*
legend
=
new
TLegend
(0.85,0.74,0.95,0.95);
legend
->SetTextFont(72);
legend
->SetTextSize(0.05);
legend
->AddEntry(
gr98
,
" q98"
,
"lp"
);
legend
->AddEntry(
gr90
,
" q90"
,
"lp"
);
legend
->AddEntry(
gr70
,
" q70"
,
"lp"
);
legend
->Draw();
}
h
#define h(i)
Definition
RSha256.hxx:106
Int_t
int Int_t
Definition
RtypesCore.h:45
Float_t
float Float_t
Definition
RtypesCore.h:57
Double_t
double Double_t
Definition
RtypesCore.h:59
kRed
@ kRed
Definition
Rtypes.h:66
kMagenta
@ kMagenta
Definition
Rtypes.h:66
kBlue
@ kBlue
Definition
Rtypes.h:66
TRangeDynCast
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Definition
TCollection.h:358
lq
int * lq
Definition
THbookFile.cxx:88
kCannotPick
@ kCannotPick
Definition
TObject.h:372
gPad
#define gPad
Definition
TVirtualPad.h:308
ROOT::Detail::TRangeCast
Definition
TCollection.h:311
TAttMarker::SetMarkerColor
virtual void SetMarkerColor(Color_t mcolor=1)
Set the marker color.
Definition
TAttMarker.h:38
TAttMarker::GetMarkerColor
virtual Color_t GetMarkerColor() const
Return the marker color.
Definition
TAttMarker.h:31
TAttMarker::SetMarkerStyle
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
Definition
TAttMarker.h:40
TAttMarker::SetMarkerSize
virtual void SetMarkerSize(Size_t msize=1)
Set the marker size.
Definition
TAttMarker.h:45
TAxis::GetXmax
Double_t GetXmax() const
Definition
TAxis.h:140
TAxis::GetXmin
Double_t GetXmin() const
Definition
TAxis.h:139
TCanvas
The Canvas class.
Definition
TCanvas.h:23
TGraph
A TGraph is an object made of two arrays X and Y with npoints each.
Definition
TGraph.h:41
TGraph::SetHighlight
virtual void SetHighlight(Bool_t set=kTRUE)
Set highlight (enable/disable) mode for the graph by default highlight mode is disable.
Definition
TGraph.cxx:2311
TGraph::Draw
void Draw(Option_t *chopt="") override
Draw this graph with its current attributes.
Definition
TGraph.cxx:831
TGraph::GetYaxis
TAxis * GetYaxis() const
Get y axis of the graph.
Definition
TGraph.cxx:1575
TGraph::SetTitle
void SetTitle(const char *title="") override
Change (i.e.
Definition
TGraph.cxx:2397
TH1F
1-D histogram with a float per channel (see TH1 documentation)
Definition
TH1.h:622
TLegend
This class displays a legend box (TPaveText) containing several legend entries.
Definition
TLegend.h:23
TList
A doubly linked list.
Definition
TList.h:38
TObject
Mother of all ROOT objects.
Definition
TObject.h:41
TString::Format
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
TText
Base class for several text objects.
Definition
TText.h:22
TVirtualPad
TVirtualPad is an abstract base class for the Pad and Canvas classes.
Definition
TVirtualPad.h:51
y
Double_t y[n]
Definition
legend1.C:17
c1
return c1
Definition
legend1.C:41
gr
TGraphErrors * gr
Definition
legend1.C:25
Authors
Rene Brun, Eddy Offermann, Jan Musinsky
Definition in file
hlquantiles.C
.
tutorials
math
hlquantiles.C
ROOT tags/6-34-04 - Reference Guide Generated on Sun Mar 23 2025 04:41:54 (GVA Time) using Doxygen 1.10.0