Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
HypoTestInverterPlot.cxx
Go to the documentation of this file.
1// @(#)root/roostats:$Id$
2// Author: Kyle Cranmer, Lorenzo Moneta, Gregory Schott, Wouter Verkerke
3/*************************************************************************
4 * Copyright (C) 1995-2008, Rene Brun and Fons Rademakers. *
5 * All rights reserved. *
6 * *
7 * For the licensing terms see $ROOTSYS/LICENSE. *
8 * For the list of contributors see $ROOTSYS/README/CREDITS. *
9 *************************************************************************/
10
11/** \class RooStats::HypoTestInverterPlot
12 \ingroup Roostats
13
14Class to plot a HypoTestInverterResult, the output of the HypoTestInverter calculator.
15
16It can be used to plot the obtained p-values ( CLb, CLs+b or CLs) for each scanned point, as well as
17the test statistic distributions (when a calculator based on pseudo-experiments is used) for the two
18hypotheses.
19
20*/
21
22// include header file of this class
24
25// include other header files
30
31#include "TGraphErrors.h"
32#include "TGraphAsymmErrors.h"
33#include "TMultiGraph.h"
34#include "TROOT.h"
35#include "TLine.h"
36#include "TAxis.h"
37#include "TLegend.h"
38#include "TVirtualPad.h"
39#include "TError.h"
41
42#include <cmath>
43
45
46using namespace RooStats;
47
48////////////////////////////////////////////////////////////////////////////////
49/// constructor from a HypoTestInverterResult class
50/// name and title are taken from the result class
51
53 TNamed( results->GetName(), results->GetTitle() ),
54 fResults(results)
55{
56}
57
58////////////////////////////////////////////////////////////////////////////////
59/// constructor with name and title from a HypoTestInverterResult class
60
62 const char* title,
63 HypoTestInverterResult* results ) :
64 TNamed( TString(name), TString(title) ),
65 fResults(results)
66{
67}
68
69////////////////////////////////////////////////////////////////////////////////
70/// Make the plot of the result of the scan using the observed data.
71/// By default plot CLs or CLsb depending if the flag UseCLs is set for the results
72/// that are passed to this instance.
73///
74/// \param opt Options according to following list:
75/// - Empty: Return CLs or CLs+b depending on the value of UseCLs.ƒ
76/// - "CLB": return CLb plot
77/// - "CLS+B" / "CLSPLUSB": return CLs+b plot independently of the flag
78/// - "CLS": return CLs plot independently of the flag
79
81{
82 TString option(opt);
83 option.ToUpper();
84 enum Mode_t { Default, CLb, CLsPlusb, CLs };
85 Mode_t type = Default;
86 if (option.Contains("CLB")) type = CLb;
87 else if (option.Contains("CLS+B") || option.Contains("CLSPLUSB")) type = CLsPlusb;
88 else if (option.Contains("CLS" )) type = CLs;
89
90 const int nEntries = fResults->ArraySize();
91
92 // sort the arrays based on the x values
93 std::vector<unsigned int> index(nEntries);
94 TMath::SortItr(fResults->fXValues.begin(), fResults->fXValues.end(), index.begin(), false);
95
96 // copy result in sorted arrays
97 std::vector<double> xArray;
98 std::vector<double> yArray;
99 std::vector<double> yErrArray;
100
101 for (int i=0; i<nEntries; i++) {
102 double CLVal = 0.;
103 double CLErr = 0.;
104 if (type == Default) {
105 CLVal = fResults->GetYValue(index[i]);
106 CLErr = fResults->GetYError(index[i]);
107 } else if (type == CLb) {
108 CLVal = fResults->CLb(index[i]);
109 CLErr = fResults->CLbError(index[i]);
110 } else if (type == CLsPlusb) {
111 CLVal = fResults->CLsplusb(index[i]);
112 CLErr = fResults->CLsplusbError(index[i]);
113 } else if (type == CLs) {
114 CLVal = fResults->CLs(index[i]);
115 CLErr = fResults->CLsError(index[i]);
116 }
117
118 if (CLVal < 0. || !std::isfinite(CLVal)) {
119 Warning("HypoTestInverterPlot::MakePlot", "Got a confidence level of %f at x=%f (failed fit?). Skipping this point.", CLVal, fResults->GetXValue(index[i]));
120 continue;
121 }
122
123 yArray.push_back(CLVal);
124 yErrArray.push_back(CLErr);
125 xArray.push_back(fResults->GetXValue(index[i]));
126 }
127
128 TGraphErrors* graph = new TGraphErrors(static_cast<Int_t>(xArray.size()),&xArray.front(),&yArray.front(),nullptr,&yErrArray.front());
129 TString pValueName = "CLs";
130 if (type == CLb ) pValueName = "CLb";
131 if (type == CLsPlusb || (type == Default && !fResults->fUseCLs) ) pValueName = "CLs+b";
132 TString name = pValueName + TString("_observed");
133 TString title = TString("Observed ") + pValueName;
134 graph->SetName(name);
135 graph->SetTitle(title);
136 graph->SetMarkerStyle(20);
137 graph->SetLineWidth(2);
138 return graph;
139}
140
141////////////////////////////////////////////////////////////////////////////////
142/// Make the expected plot and the bands
143/// nsig1 and nsig2 indicates the n-sigma value for the bands
144/// if nsig1 = 0 no band is drawn (only expected value)
145/// if nsig2 > nsig1 (default is nsig1=1 and nsig2=2) the second band is also drawn
146/// The first band is drawn in green while the second in yellow
147/// THe return result is a TMultiGraph object
148
150{
151 const int nEntries = fResults->ArraySize();
152 bool doFirstBand = (nsig1 > 0);
153 bool doSecondBand = (nsig2 > nsig1);
154
155 nsig1 = std::abs(nsig1);
156 nsig2 = std::abs(nsig2);
157
158 // sort the arrays based on the x values
159 std::vector<unsigned int> index(nEntries);
160 TMath::SortItr(fResults->fXValues.begin(), fResults->fXValues.end(), index.begin(), false);
161
162 // create the graphs
163 TGraph * g0 = new TGraph;
164 TString pValueName = "CLs";
165 if (!fResults->fUseCLs) pValueName = "CLs+b";
166 g0->SetTitle(TString::Format("Expected %s - Median",pValueName.Data()) );
167 TGraphAsymmErrors * g1 = nullptr;
168 TGraphAsymmErrors * g2 = nullptr;
169 if (doFirstBand) {
170 g1 = new TGraphAsymmErrors;
171 if (nsig1 - int(nsig1) < 0.01) {
172 g1->SetTitle(TString::Format("Expected %s #pm %d #sigma",pValueName.Data(),int(nsig1)) );
173 } else {
174 g1->SetTitle(TString::Format("Expected %s #pm %3.1f #sigma", pValueName.Data(), nsig1));
175 }
176 }
177 if (doSecondBand) {
178 g2 = new TGraphAsymmErrors;
179 if (nsig2 - int(nsig2) < 0.01) {
180 g2->SetTitle(TString::Format("Expected %s #pm %d #sigma",pValueName.Data(),int(nsig2)) );
181 } else {
182 g2->SetTitle(TString::Format("Expected %s #pm %3.1f #sigma", pValueName.Data(), nsig2));
183 }
184 }
185 double p[7];
186 double q[7];
187 p[0] = ROOT::Math::normal_cdf(-nsig2);
188 p[1] = ROOT::Math::normal_cdf(-nsig1);
189 p[2] = 0.5;
190 p[3] = ROOT::Math::normal_cdf(nsig1);
191 p[4] = ROOT::Math::normal_cdf(nsig2);
192
193 bool resultIsAsymptotic = ( !fResults->GetNullTestStatDist(0) && !fResults->GetAltTestStatDist(0) );
194 int np = 0;
195 for (int j=0; j<nEntries; ++j) {
196 int i = index[j]; // i is the order index
198 if ( !s) continue;
199 const std::vector<double> & values = s->GetSamplingDistribution();
200 // special case for asymptotic results (cannot use TMath::quantile in that case)
201 if (resultIsAsymptotic) {
202 double maxSigma = fResults->fgAsymptoticMaxSigma;
203 double dsig = 2* maxSigma/ (values.size() -1) ;
204 int i0 = (int) TMath::Floor ( ( -nsig2 + maxSigma )/dsig + 0.5);
205 int i1 = (int) TMath::Floor ( (-nsig1 + maxSigma )/dsig + 0.5);
206 int i2 = (int) TMath::Floor ( ( maxSigma)/dsig + 0.5);
207 int i3 = (int) TMath::Floor ( ( nsig1 + maxSigma)/dsig + 0.5);
208 int i4 = (int) TMath::Floor ( ( nsig2 + maxSigma)/dsig + 0.5);
209 q[0] = values[i0];
210 q[1] = values[i1];
211 q[2] = values[i2];
212 q[3] = values[i3];
213 q[4] = values[i4];
214 }
215 else {
216 double * x = const_cast<double *>(&values[0]); // need to change TMath::Quantiles
217 TMath::Quantiles(values.size(), 5, x,q,p,false);
218 }
219
220 g0->SetPoint(np, fResults->GetXValue(i), q[2]);
221 if (g1) {
222 g1->SetPoint(np, fResults->GetXValue(i), q[2]);
223 g1->SetPointEYlow(np, q[2] - q[1]); // -1 sigma error
224 g1->SetPointEYhigh(np, q[3] - q[2]);//+1 sigma error
225 }
226 if (g2) {
227 g2->SetPoint(np, fResults->GetXValue(i), q[2]);
228
229 g2->SetPointEYlow(np, q[2]-q[0]); // -2 -- -1 sigma error
230 g2->SetPointEYhigh(np, q[4]-q[2]);
231 }
232 if (s) delete s;
233 np++;
234 }
235
236 TString name = GetName() + TString("_expected");
237 TString title = TString("Expected ") + GetTitle();
238 TMultiGraph* graph = new TMultiGraph(name,title);
239
240 // set the graphics options and add in multi graph
241 // orderof adding is drawing order
242 if (g2) {
244 graph->Add(g2,"3");
245 }
246 if (g1) {
247 g1->SetFillColor(kGreen);
248 graph->Add(g1,"3");
249 }
250 g0->SetLineStyle(2);
251 g0->SetLineWidth(2);
252 graph->Add(g0,"L");
253
254 return graph;
255}
256
257////////////////////////////////////////////////////////////////////////////////
258/// destructor
259
261{
262}
263
264////////////////////////////////////////////////////////////////////////////////
265/// Draw the result in the current canvas
266/// Possible options:
267/// - SAME : draw in the current axis
268/// - OBS : draw only the observed plot
269/// - EXP : draw only the expected plot
270///
271/// - CLB : draw also the CLB
272/// - 2CL : draw both clsplusb and cls
273///
274/// default draw observed + expected with 1 and 2 sigma bands
275
277
278 TString option(opt);
279 option.ToUpper();
280 bool drawAxis = !option.Contains("SAME");
281 bool drawObs = option.Contains("OBS") || !option.Contains("EXP");
282 bool drawExp = option.Contains("EXP") || !option.Contains("OBS");
283 bool drawCLb = option.Contains("CLB");
284 bool draw2CL = option.Contains("2CL");
285
286 TGraphErrors * gobs = nullptr;
287 TGraph * gplot = nullptr;
288 if (drawObs) {
289 gobs = MakePlot();
290 // add object to top-level directory to avoid mem leak
291 if (gROOT) gROOT->Add(gobs);
292 if (drawAxis) {
293 gobs->Draw("APL");
294 gplot = gobs;
295 gplot->GetHistogram()->SetTitle( GetTitle() );
296 }
297 else gobs->Draw("PL");
298
299 }
300 TMultiGraph * gexp = nullptr;
301 if (drawExp) {
302 gexp = MakeExpectedPlot();
303 // add object to current directory to avoid mem leak
304 if (gROOT) gROOT->Add(gexp);
305 if (drawAxis && !drawObs) {
306 gexp->Draw("A");
307 if (gexp->GetHistogram()) gexp->GetHistogram()->SetTitle( GetTitle() );
308 gplot = static_cast<TGraph*>(gexp->GetListOfGraphs()->First());
309 }
310 else
311 gexp->Draw();
312
313 }
314
315 // draw also an horizontal line at the desired conf level
316 if (gplot) {
317 double alpha = 1.-fResults->ConfidenceLevel();
318 double x1 = gplot->GetXaxis()->GetXmin();
319 double x2 = gplot->GetXaxis()->GetXmax();
320 TLine * line = new TLine(x1, alpha, x2,alpha);
322 line->Draw();
323 // put axis labels
325 if (arg) gplot->GetXaxis()->SetTitle(arg->GetName());
326 gplot->GetYaxis()->SetTitle("p value");
327 }
328
329
330 TGraph *gclb = nullptr;
331 if (drawCLb) {
332 gclb = MakePlot("CLb");
333 if (gROOT) gROOT->Add(gclb);
334 gclb->SetMarkerColor(kBlue+4);
335 gclb->Draw("PL");
336 // draw in red observed cls or clsb
337 if (gobs) gobs->SetMarkerColor(kRed);
338 }
339 TGraph * gclsb = nullptr;
340 TGraph * gcls = nullptr;
341 if (draw2CL) {
342 if (fResults->fUseCLs) {
343 gclsb = MakePlot("CLs+b");
344 if (gROOT) gROOT->Add(gclsb);
345 gclsb->SetMarkerColor(kBlue);
346 gclsb->Draw("PL");
347 gclsb->SetLineStyle(3);
348 }
349 else {
350 gcls = MakePlot("CLs");
351 if (gROOT) gROOT->Add(gcls);
352 gcls->SetMarkerColor(kBlue);
353 gcls->Draw("PL");
354 gcls->SetLineStyle(3);
355 }
356 }
357 // draw again observed values otherwise will be covered by the bands
358 if (gobs) {
359 gobs->Draw("PL");
360 }
361
362
363 double y0 = 0.6;
364 double verticalSize = (gexp || draw2CL || drawCLb ) ? 0.3 : 0.15;
365 double y1 = y0 + verticalSize;
366 TLegend * l = new TLegend(0.6,y0,0.9,y1);
367 if (gobs) l->AddEntry(gobs,"","PEL");
368 if (gclsb) l->AddEntry(gclsb,"","PEL");
369 if (gcls) l->AddEntry(gcls,"","PEL");
370 if (gclb) l->AddEntry(gclb,"","PEL");
371 if (gexp) {
372 // loop in reverse order (opposite to drawing one)
373 int ngraphs = gexp->GetListOfGraphs()->GetSize();
374 for (int i = ngraphs-1; i>=0; --i) {
375 TObject * obj = gexp->GetListOfGraphs()->At(i);
376 TString lopt = "F";
377 if (i == ngraphs-1) lopt = "L";
378 if (obj) l->AddEntry(obj,"",lopt);
379 }
380 }
381 l->Draw();
382 // redraw the axis
383 if (gPad) gPad->RedrawAxis();
384
385}
386
387////////////////////////////////////////////////////////////////////////////////
388/// \param index Index of the result stored in HypoTestInverterResult
389/// \param type Type of the test (see below)
390/// \param nbins Number of bins
391/// - type =0 null and alt
392/// - type = 1 only null (S+B)
393/// - type = 2 only alt (B)
394
396 SamplingDistPlot * pl = nullptr;
397 if (type == 0) {
399 if (result)
400 pl = new HypoTestPlot(*result, nbins );
401 return pl;
402 }
403 if (type == 1) {
405 if (sbDist) {
406 pl = new SamplingDistPlot( nbins);
407 pl->AddSamplingDistribution(sbDist);
408 return pl;
409 }
410 }
411 if (type == 2) {
413 if (bDist) {
414 pl = new SamplingDistPlot( nbins);
415 pl->AddSamplingDistribution(bDist);
416 return pl;
417 }
418 }
419 return nullptr;
420}
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:377
@ kRed
Definition Rtypes.h:66
@ kGreen
Definition Rtypes.h:66
@ kBlue
Definition Rtypes.h:66
@ kYellow
Definition Rtypes.h:66
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t option
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t np
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t result
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t index
Option_t Option_t TPoint TPoint const char x2
Option_t Option_t TPoint TPoint const char x1
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t type
Option_t Option_t TPoint TPoint const char y1
char name[80]
Definition TGX11.cxx:110
float * q
#define gROOT
Definition TROOT.h:406
#define gPad
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:77
RooAbsArg * first() const
Class to plot a HypoTestInverterResult, the output of the HypoTestInverter calculator.
HypoTestInverterResult * fResults
void Draw(Option_t *opt="") override
Draw the scan result in the current canvas Possible options: "" (default): draw observed + expected w...
SamplingDistPlot * MakeTestStatPlot(int index, int type=0, int nbins=100)
Plot the test statistic distributions.
HypoTestInverterPlot(HypoTestInverterResult *results)
constructor
TMultiGraph * MakeExpectedPlot(double sig1=1, double sig2=2)
Make the expected plot and the bands nsig1 and nsig2 indicates the n-sigma value for the bands if nsi...
TGraphErrors * MakePlot(Option_t *opt="")
return a TGraphErrors with the obtained observed p-values resultinf from the scan By default (Option ...
HypoTestInverterResult class holds the array of hypothesis test results and compute a confidence inte...
double GetYValue(int index) const
function to return the value of the confidence level for the i^th entry in the results
SamplingDistribution * GetSignalAndBackgroundTestStatDist(int index) const
get the signal and background test statistic distribution
int ArraySize() const
number of entries in the results array
double CLsplusbError(int index) const
return the observed CLsplusb value for the i-th entry
double CLsError(int index) const
return the observed CLb value for the i-th entry
double CLbError(int index) const
return the observed CLb value for the i-th entry
double GetXValue(int index) const
function to return the value of the parameter of interest for the i^th entry in the results
static double fgAsymptoticMaxSigma
max sigma value used to scan asymptotic expected p values
double CLs(int index) const
return the observed CLb value for the i-th entry
double GetYError(int index) const
function to return the estimated error on the value of the confidence level for the i^th entry in the...
double CLb(int index) const
return the observed CLb value for the i-th entry
TList fYObjects
list of HypoTestResult for each point
SamplingDistribution * GetExpectedPValueDist(int index) const
return expected distribution of p-values (Cls or Clsplusb)
double CLsplusb(int index) const
return the observed CLsplusb value for the i-th entry
SamplingDistribution * GetNullTestStatDist(int index) const
same in terms of alt and null
SamplingDistribution * GetAltTestStatDist(int index) const
SamplingDistribution * GetBackgroundTestStatDist(int index) const
get the background test statistic distribution
This class provides the plots for the result of a study performed with any of the HypoTestCalculatorG...
HypoTestResult is a base class for results from hypothesis tests.
This class provides simple and straightforward utilities to plot SamplingDistribution objects.
double AddSamplingDistribution(const SamplingDistribution *samplingDist, Option_t *drawOptions="NORMALIZE HIST")
adds the sampling distribution and returns the scale factor
This class simply holds a sampling distribution of some test statistic.
const std::vector< double > & GetSamplingDistribution() const
Get test statistics values.
RooArgSet fParameters
set containing the parameter of interest
double ConfidenceLevel() const override
return the confidence interval
virtual void SetFillColor(Color_t fcolor)
Set the fill area color.
Definition TAttFill.h:37
virtual void SetLineStyle(Style_t lstyle)
Set the line style.
Definition TAttLine.h:42
virtual void SetLineWidth(Width_t lwidth)
Set the line width.
Definition TAttLine.h:43
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition TAttLine.h:40
virtual void SetMarkerColor(Color_t mcolor=1)
Set the marker color.
Definition TAttMarker.h:38
Double_t GetXmax() const
Definition TAxis.h:140
Double_t GetXmin() const
Definition TAxis.h:139
virtual Int_t GetSize() const
Return the capacity of the collection, i.e.
TGraph with asymmetric error bars.
virtual void SetPointEYlow(Int_t i, Double_t eyl)
Set EYlow for point i.
virtual void SetPointEYhigh(Int_t i, Double_t eyh)
Set EYhigh for point i.
A TGraphErrors is a TGraph with error bars.
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
void Draw(Option_t *chopt="") override
Draw this graph with its current attributes.
Definition TGraph.cxx:809
TAxis * GetXaxis() const
Get x axis of the graph.
Definition TGraph.cxx:1544
TAxis * GetYaxis() const
Get y axis of the graph.
Definition TGraph.cxx:1553
virtual TH1F * GetHistogram() const
Returns a pointer to the histogram used to draw the axis Takes into account the two following cases.
Definition TGraph.cxx:1406
void SetTitle(const char *title="") override
Change (i.e.
Definition TGraph.cxx:2374
void SetTitle(const char *title) override
Change/set the title.
Definition TH1.cxx:6686
This class displays a legend box (TPaveText) containing several legend entries.
Definition TLegend.h:23
Use the TLine constructor to create a simple line.
Definition TLine.h:22
TObject * First() const override
Return the first object in the list. Returns 0 when list is empty.
Definition TList.cxx:657
TObject * At(Int_t idx) const override
Returns the object at position idx. Returns 0 if idx is out of range.
Definition TList.cxx:355
A TMultiGraph is a collection of TGraph (or derived) objects.
Definition TMultiGraph.h:34
TList * GetListOfGraphs() const
Definition TMultiGraph.h:68
TH1F * GetHistogram()
Returns a pointer to the histogram used to draw the axis.
void Draw(Option_t *chopt="") override
Draw this multigraph with its current attributes.
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition TNamed.cxx:164
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
const char * GetTitle() const override
Returns title of object.
Definition TNamed.h:48
Mother of all ROOT objects.
Definition TObject.h:41
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:973
virtual void Draw(Option_t *option="")
Default Draw method for all objects.
Definition TObject.cxx:274
Basic string class.
Definition TString.h:139
const char * Data() const
Definition TString.h:376
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
TLine * line
double normal_cdf(double x, double sigma=1, double x0=0)
Cumulative distribution function of the normal (Gaussian) distribution (lower tail).
Double_t x[n]
Definition legend1.C:17
Namespace for the RooStats classes.
Definition Asimov.h:19
void SortItr(Iterator first, Iterator last, IndexIterator index, Bool_t down=kTRUE)
Sort the n1 elements of the Short_t array defined by its iterators.
Definition TMathBase.h:406
Double_t Floor(Double_t x)
Rounds x downward, returning the largest integral value that is not greater than x.
Definition TMath.h:680
void Quantiles(Int_t n, Int_t nprob, Double_t *x, Double_t *quantiles, Double_t *prob, Bool_t isSorted=kTRUE, Int_t *index=nullptr, Int_t type=7)
Computes sample quantiles, corresponding to the given probabilities.
Definition TMath.cxx:1207
Definition graph.py:1
TLine l
Definition textangle.C:4