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
29
30#include "TGraphErrors.h"
31#include "TGraphAsymmErrors.h"
32#include "TMultiGraph.h"
33#include "TROOT.h"
34#include "TLine.h"
35#include "TAxis.h"
36#include "TLegend.h"
37#include "TVirtualPad.h"
38#include "TError.h"
40
41#include <cmath>
42
43
44using namespace RooStats;
45
46////////////////////////////////////////////////////////////////////////////////
47/// constructor from a HypoTestInverterResult class
48/// name and title are taken from the result class
49
55
56////////////////////////////////////////////////////////////////////////////////
57/// constructor with name and title from a HypoTestInverterResult class
58
60 const char* title,
62 TNamed( TString(name), TString(title) ),
63 fResults(results)
64{
65}
66
67////////////////////////////////////////////////////////////////////////////////
68/// Make the plot of the result of the scan using the observed data.
69/// By default plot CLs or CLsb depending if the flag UseCLs is set for the results
70/// that are passed to this instance.
71///
72/// \param opt Options according to following list:
73/// - Empty: Return CLs or CLs+b depending on the value of UseCLs.ƒ
74/// - "CLB": return CLb plot
75/// - "CLS+B" / "CLSPLUSB": return CLs+b plot independently of the flag
76/// - "CLS": return CLs plot independently of the flag
77
79{
80 TString option(opt);
81 option.ToUpper();
82 enum Mode_t { Default, CLb, CLsPlusb, CLs };
83 Mode_t type = Default;
84 if (option.Contains("CLB")) type = CLb;
85 else if (option.Contains("CLS+B") || option.Contains("CLSPLUSB")) type = CLsPlusb;
86 else if (option.Contains("CLS" )) type = CLs;
87
88 const int nEntries = fResults->ArraySize();
89
90 // sort the arrays based on the x values
91 std::vector<unsigned int> index(nEntries);
92 TMath::SortItr(fResults->fXValues.begin(), fResults->fXValues.end(), index.begin(), false);
93
94 // copy result in sorted arrays
95 std::vector<double> xArray;
96 std::vector<double> yArray;
97 std::vector<double> yErrArray;
98
99 for (int i=0; i<nEntries; i++) {
100 double CLVal = 0.;
101 double CLErr = 0.;
102 if (type == Default) {
105 } else if (type == CLb) {
106 CLVal = fResults->CLb(index[i]);
108 } else if (type == CLsPlusb) {
111 } else if (type == CLs) {
112 CLVal = fResults->CLs(index[i]);
114 }
115
116 if (CLVal < 0. || !std::isfinite(CLVal)) {
117 Warning("HypoTestInverterPlot::MakePlot", "Got a confidence level of %f at x=%f (failed fit?). Skipping this point.", CLVal, fResults->GetXValue(index[i]));
118 continue;
119 }
120
121 yArray.push_back(CLVal);
122 yErrArray.push_back(CLErr);
123 xArray.push_back(fResults->GetXValue(index[i]));
124 }
125
126 TGraphErrors* graph = new TGraphErrors(static_cast<Int_t>(xArray.size()),&xArray.front(),&yArray.front(),nullptr,&yErrArray.front());
127 TString pValueName = "CLs";
128 if (type == CLb ) pValueName = "CLb";
129 if (type == CLsPlusb || (type == Default && !fResults->fUseCLs) ) pValueName = "CLs+b";
130 TString name = pValueName + TString("_observed");
131 TString title = TString("Observed ") + pValueName;
132 graph->SetName(name);
133 graph->SetTitle(title);
134 graph->SetMarkerStyle(20);
135 graph->SetLineWidth(2);
136 return graph;
137}
138
139////////////////////////////////////////////////////////////////////////////////
140/// Make the expected plot and the bands
141/// nsig1 and nsig2 indicates the n-sigma value for the bands
142/// if nsig1 = 0 no band is drawn (only expected value)
143/// if nsig2 > nsig1 (default is nsig1=1 and nsig2=2) the second band is also drawn
144/// The first band is drawn in green while the second in yellow
145/// THe return result is a TMultiGraph object
146
148{
149 const int nEntries = fResults->ArraySize();
150 bool doFirstBand = (nsig1 > 0);
151 bool doSecondBand = (nsig2 > nsig1);
152
153 nsig1 = std::abs(nsig1);
154 nsig2 = std::abs(nsig2);
155
156 // sort the arrays based on the x values
157 std::vector<unsigned int> index(nEntries);
158 TMath::SortItr(fResults->fXValues.begin(), fResults->fXValues.end(), index.begin(), false);
159
160 // create the graphs
161 TGraph * g0 = new TGraph;
162 TString pValueName = "CLs";
163 if (!fResults->fUseCLs) pValueName = "CLs+b";
164 g0->SetTitle(TString::Format("Expected %s - Median",pValueName.Data()) );
165 TGraphAsymmErrors * g1 = nullptr;
166 TGraphAsymmErrors * g2 = nullptr;
167 if (doFirstBand) {
168 g1 = new TGraphAsymmErrors;
169 if (nsig1 - int(nsig1) < 0.01) {
170 g1->SetTitle(TString::Format("Expected %s #pm %d #sigma",pValueName.Data(),int(nsig1)) );
171 } else {
172 g1->SetTitle(TString::Format("Expected %s #pm %3.1f #sigma", pValueName.Data(), nsig1));
173 }
174 }
175 if (doSecondBand) {
176 g2 = new TGraphAsymmErrors;
177 if (nsig2 - int(nsig2) < 0.01) {
178 g2->SetTitle(TString::Format("Expected %s #pm %d #sigma",pValueName.Data(),int(nsig2)) );
179 } else {
180 g2->SetTitle(TString::Format("Expected %s #pm %3.1f #sigma", pValueName.Data(), nsig2));
181 }
182 }
183 double p[7];
184 double q[7];
187 p[2] = 0.5;
190
192 int np = 0;
193 for (int j=0; j<nEntries; ++j) {
194 int i = index[j]; // i is the order index
196 if ( !s) continue;
197 const std::vector<double> & values = s->GetSamplingDistribution();
198 // special case for asymptotic results (cannot use TMath::quantile in that case)
199 if (resultIsAsymptotic) {
201 double dsig = 2* maxSigma/ (values.size() -1) ;
202 int i0 = (int) TMath::Floor ( ( -nsig2 + maxSigma )/dsig + 0.5);
203 int i1 = (int) TMath::Floor ( (-nsig1 + maxSigma )/dsig + 0.5);
204 int i2 = (int) TMath::Floor ( ( maxSigma)/dsig + 0.5);
205 int i3 = (int) TMath::Floor ( ( nsig1 + maxSigma)/dsig + 0.5);
206 int i4 = (int) TMath::Floor ( ( nsig2 + maxSigma)/dsig + 0.5);
207 q[0] = values[i0];
208 q[1] = values[i1];
209 q[2] = values[i2];
210 q[3] = values[i3];
211 q[4] = values[i4];
212 }
213 else {
214 double * x = const_cast<double *>(&values[0]); // need to change TMath::Quantiles
215 TMath::Quantiles(values.size(), 5, x,q,p,false);
216 }
217
218 g0->SetPoint(np, fResults->GetXValue(i), q[2]);
219 if (g1) {
220 g1->SetPoint(np, fResults->GetXValue(i), q[2]);
221 g1->SetPointEYlow(np, q[2] - q[1]); // -1 sigma error
222 g1->SetPointEYhigh(np, q[3] - q[2]);//+1 sigma error
223 }
224 if (g2) {
225 g2->SetPoint(np, fResults->GetXValue(i), q[2]);
226
227 g2->SetPointEYlow(np, q[2]-q[0]); // -2 -- -1 sigma error
228 g2->SetPointEYhigh(np, q[4]-q[2]);
229 }
230 if (s) delete s;
231 np++;
232 }
233
234 TString name = GetName() + TString("_expected");
235 TString title = TString("Expected ") + GetTitle();
236 TMultiGraph* graph = new TMultiGraph(name,title);
237
238 // set the graphics options and add in multi graph
239 // orderof adding is drawing order
240 if (g2) {
241 g2->SetFillColor(kYellow);
242 graph->Add(g2,"3");
243 }
244 if (g1) {
245 g1->SetFillColor(kGreen);
246 graph->Add(g1,"3");
247 }
248 g0->SetLineStyle(2);
249 g0->SetLineWidth(2);
250 graph->Add(g0,"L");
251
252 return graph;
253}
254
255////////////////////////////////////////////////////////////////////////////////
256/// destructor
257
261
262////////////////////////////////////////////////////////////////////////////////
263/// Draw the result in the current canvas
264/// Possible options:
265/// - SAME : draw in the current axis
266/// - OBS : draw only the observed plot
267/// - EXP : draw only the expected plot
268///
269/// - CLB : draw also the CLB
270/// - 2CL : draw both clsplusb and cls
271///
272/// default draw observed + expected with 1 and 2 sigma bands
273
275
276 TString option(opt);
277 option.ToUpper();
278 bool drawAxis = !option.Contains("SAME");
279 bool drawObs = option.Contains("OBS") || !option.Contains("EXP");
280 bool drawExp = option.Contains("EXP") || !option.Contains("OBS");
281 bool drawCLb = option.Contains("CLB");
282 bool draw2CL = option.Contains("2CL");
283
284 TGraphErrors * gobs = nullptr;
285 TGraph * gplot = nullptr;
286 if (drawObs) {
287 gobs = MakePlot();
288 // add object to top-level directory to avoid mem leak
289 if (gROOT) gROOT->Add(gobs);
290 if (drawAxis) {
291 gobs->Draw("APL");
292 gplot = gobs;
293 gplot->GetHistogram()->SetTitle( GetTitle() );
294 }
295 else gobs->Draw("PL");
296
297 }
298 TMultiGraph * gexp = nullptr;
299 if (drawExp) {
301 // add object to current directory to avoid mem leak
302 if (gROOT) gROOT->Add(gexp);
303 if (drawAxis && !drawObs) {
304 gexp->Draw("A");
305 if (gexp->GetHistogram()) gexp->GetHistogram()->SetTitle( GetTitle() );
306 gplot = static_cast<TGraph*>(gexp->GetListOfGraphs()->First());
307 }
308 else
309 gexp->Draw();
310
311 }
312
313 // draw also an horizontal line at the desired conf level
314 if (gplot) {
315 double alpha = 1.-fResults->ConfidenceLevel();
316 double x1 = gplot->GetXaxis()->GetXmin();
317 double x2 = gplot->GetXaxis()->GetXmax();
318 TLine * line = new TLine(x1, alpha, x2,alpha);
320 line->Draw();
321 // put axis labels
323 if (arg) gplot->GetXaxis()->SetTitle(arg->GetName());
324 gplot->GetYaxis()->SetTitle("p value");
325 }
326
327
328 TGraph *gclb = nullptr;
329 if (drawCLb) {
330 gclb = MakePlot("CLb");
331 if (gROOT) gROOT->Add(gclb);
332 gclb->SetMarkerColor(kBlue+4);
333 gclb->Draw("PL");
334 // draw in red observed cls or clsb
335 if (gobs) gobs->SetMarkerColor(kRed);
336 }
337 TGraph * gclsb = nullptr;
338 TGraph * gcls = nullptr;
339 if (draw2CL) {
340 if (fResults->fUseCLs) {
341 gclsb = MakePlot("CLs+b");
342 if (gROOT) gROOT->Add(gclsb);
343 gclsb->SetMarkerColor(kBlue);
344 gclsb->Draw("PL");
345 gclsb->SetLineStyle(3);
346 }
347 else {
348 gcls = MakePlot("CLs");
349 if (gROOT) gROOT->Add(gcls);
350 gcls->SetMarkerColor(kBlue);
351 gcls->Draw("PL");
352 gcls->SetLineStyle(3);
353 }
354 }
355 // draw again observed values otherwise will be covered by the bands
356 if (gobs) {
357 gobs->Draw("PL");
358 }
359
360
361 double y0 = 0.6;
362 double verticalSize = (gexp || draw2CL || drawCLb ) ? 0.3 : 0.15;
363 double y1 = y0 + verticalSize;
364 TLegend * l = new TLegend(0.6,y0,0.9,y1);
365 if (gobs) l->AddEntry(gobs,"","PEL");
366 if (gclsb) l->AddEntry(gclsb,"","PEL");
367 if (gcls) l->AddEntry(gcls,"","PEL");
368 if (gclb) l->AddEntry(gclb,"","PEL");
369 if (gexp) {
370 // loop in reverse order (opposite to drawing one)
371 int ngraphs = gexp->GetListOfGraphs()->GetSize();
372 for (int i = ngraphs-1; i>=0; --i) {
373 TObject * obj = gexp->GetListOfGraphs()->At(i);
374 TString lopt = "F";
375 if (i == ngraphs-1) lopt = "L";
376 if (obj) l->AddEntry(obj,"",lopt);
377 }
378 }
379 l->Draw();
380 // redraw the axis
381 if (gPad) gPad->RedrawAxis();
382
383}
384
385////////////////////////////////////////////////////////////////////////////////
386/// \param index Index of the result stored in HypoTestInverterResult
387/// \param type Type of the test (see below)
388/// \param nbins Number of bins
389/// - type =0 null and alt
390/// - type = 1 only null (S+B)
391/// - type = 2 only alt (B)
392
394 SamplingDistPlot * pl = nullptr;
395 if (type == 0) {
397 if (result)
398 pl = new HypoTestPlot(*result, nbins );
399 return pl;
400 }
401 if (type == 1) {
403 if (sbDist) {
404 pl = new SamplingDistPlot( nbins);
405 pl->AddSamplingDistribution(sbDist);
406 return pl;
407 }
408 }
409 if (type == 2) {
411 if (bDist) {
412 pl = new SamplingDistPlot( nbins);
413 pl->AddSamplingDistribution(bDist);
414 return pl;
415 }
416 }
417 return nullptr;
418}
const char Option_t
Option string (const char)
Definition RtypesCore.h:80
@ kRed
Definition Rtypes.h:67
@ kGreen
Definition Rtypes.h:67
@ kBlue
Definition Rtypes.h:67
@ kYellow
Definition Rtypes.h:67
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
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:146
float * q
#define gROOT
Definition TROOT.h:414
#define gPad
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:76
RooAbsArg * first() const
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.
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 SetLineWidth(Width_t lwidth)
Set the line width.
Definition TAttLine.h:47
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition TAttLine.h:44
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
Definition TAttMarker.h:43
TGraph with asymmetric error bars.
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
void SetName(const char *name="") override
Set graph name.
Definition TGraph.cxx:2426
void SetTitle(const char *title="") override
Change (i.e.
Definition TGraph.cxx:2442
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 * At(Int_t idx) const override
Returns the object at position idx. Returns 0 if idx is out of range.
Definition TList.cxx:487
A TMultiGraph is a collection of TGraph (or derived) objects.
Definition TMultiGraph.h:34
virtual void Add(TGraph *graph, Option_t *chopt="")
Add a new graph to the list of graphs.
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
const char * GetTitle() const override
Returns title of object.
Definition TNamed.h:50
Mother of all ROOT objects.
Definition TObject.h:42
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:1081
virtual void Draw(Option_t *option="")
Default Draw method for all objects.
Definition TObject.cxx:290
Basic string class.
Definition TString.h:138
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:2384
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 CodegenImpl.h:65
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:388
Double_t Floor(Double_t x)
Rounds x downward, returning the largest integral value that is not greater than x.
Definition TMath.h:691
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
TLine l
Definition textangle.C:4