Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
Loading...
Searching...
No Matches
HybridPlot.cxx
Go to the documentation of this file.
1// @(#)root/hist:$Id$
2
3/** \class RooStats::HybridPlot
4 \ingroup Roostats
5
6This class provides the plots for the result of a study performed with the
7HybridCalculatorOriginal class.
8
9Authors: D. Piparo, G. Schott - Universitaet Karlsruhe
10
11An example plot is available here:
12http://www-ekp.physik.uni-karlsruhe.de/~schott/roostats/hybridplot_example.png
13*/
14
15#include <cassert>
16#include <cmath>
17#include <iostream>
18#include <map>
19
20#include "RooStats/HybridPlot.h"
21#include "TStyle.h"
22#include "TF1.h"
23#include "TAxis.h"
24#include "TH1.h"
25#include "TLine.h"
26#include "TLegend.h"
27#include "TFile.h"
28#include "TVirtualPad.h"
29
30#include <algorithm>
31
32/// To build the THtml documentation
33
34
35using namespace RooStats;
36
37////////////////////////////////////////////////////////////////////////////////
38/// HybridPlot constructor
39
41 const char* title,
42 const std::vector<double> & sb_vals,
43 const std::vector<double> & b_vals,
44 double testStat_data,
45 int n_bins,
46 bool verbosity):
47 TNamed(name,title),
48 fSb_histo(nullptr),
49 fSb_histo_shaded(nullptr),
50 fB_histo(nullptr),
51 fB_histo_shaded(nullptr),
52 fData_testStat_line(nullptr),
53 fLegend(nullptr),
54 fPad(nullptr),
56{
57 int nToysSB = sb_vals.size();
58 int nToysB = sb_vals.size();
59 assert (nToysSB >0);
60 assert (nToysB >0);
61
62 // Get the max and the min of the plots
63 double min = *std::min_element(sb_vals.begin(), sb_vals.end());
64 double max = *std::max_element(sb_vals.begin(), sb_vals.end());
65
66 double min_b = *std::min_element(b_vals.begin(), b_vals.end());
67 double max_b = *std::max_element(b_vals.begin(), b_vals.end());
68
69
70 if ( min_b < min) min = min_b;
71 if ( max_b > max) max = max_b;
72
73 if (testStat_data<min) min = testStat_data;
74 if (testStat_data>max) max = testStat_data;
75
76 min *= 1.1;
77 max *= 1.1;
78
79 // Build the histos
80
81 fSb_histo = new TH1F ("SB_model",title,n_bins,min,max);
84 fSb_histo->GetXaxis()->SetTitle("test statistics");
86
87 fB_histo = new TH1F ("B_model",title,n_bins,min,max);
90 fB_histo->GetXaxis()->SetTitle("test statistics");
92
93 for (int i=0;i<nToysSB;++i) fSb_histo->Fill(sb_vals[i]);
94 for (int i=0;i<nToysB;++i) fB_histo->Fill(b_vals[i]);
95
99
100 // Build the line of the measured -2lnQ
104
105 // The legend
106 double golden_section = (std::sqrt(5.)-1)/2;
107
108 fLegend = new TLegend(0.75,0.95-0.2*golden_section,0.95,0.95);
109 TString title_leg="test statistics distributions ";
110 title_leg+=sb_vals.size();
111 title_leg+=" toys";
112 fLegend->SetName(title_leg.Data());
113 fLegend->AddEntry(fSb_histo,"SB toy datasets");
114 fLegend->AddEntry(fB_histo,"B toy datasets");
115 fLegend->AddEntry((TLine*)fData_testStat_line,"test statistics on data","L");
117}
118
119////////////////////////////////////////////////////////////////////////////////
120/// destructor
121
123{
124
125 if (fSb_histo) delete fSb_histo;
126 if (fB_histo) delete fB_histo;
130 if (fLegend) delete fLegend;
131}
132
133////////////////////////////////////////////////////////////////////////////////
134/// draw the S+B and B histogram in the current canvas
135
136void HybridPlot::Draw(const char* )
137{
138 // We don't want the statistics of the histos
139 gStyle->SetOptStat(0);
140
141 // The histos
142
145 fB_histo->DrawNormalized("same");
146 }
147 else{
149 fSb_histo->DrawNormalized("same");
150 }
151
152 // Shaded
153 fB_histo_shaded = static_cast<TH1F*>(fB_histo->Clone("b_shaded"));
156
157 fSb_histo_shaded = static_cast<TH1F*>(fSb_histo->Clone("sb_shaded"));
160
161 // Empty the bins according to the data -2lnQ
163 for (int i=1;i<=fSb_histo->GetNbinsX();++i){
167 }
168 else{
171 }
172 }
173
174 // Draw the shaded histos
175 fB_histo_shaded->Draw("same");
176 fSb_histo_shaded->Draw("same");
177
178 // The line
179 fData_testStat_line->Draw("same");
180
181 // The legend
182 fLegend->Draw("same");
183
184 if (gPad) {
185 gPad->SetName(GetName());
186 gPad->SetTitle(GetTitle() );
187 }
188
189 fPad = gPad;
190
191}
192
193////////////////////////////////////////////////////////////////////////////////
194/// All the objects are written to rootfile
195
196void HybridPlot::DumpToFile (const char* RootFileName, const char* options)
197{
198
199 TFile ofile(RootFileName,options);
200 ofile.cd();
201
202 // The histos
203 fSb_histo->Write();
204 fB_histo->Write();
205
206 // The shaded histos
207 if (fB_histo_shaded!=nullptr && fSb_histo_shaded!=nullptr){
210 }
211
212 // The line
213 fData_testStat_line->Write("Measured test statistics line tag");
214
215 // The legend
216 fLegend->Write();
217
218 ofile.Close();
219
220}
221
222////////////////////////////////////////////////////////////////////////////////
223
225 if (!fPad) {
226 Error("HybridPlot","Hybrid plot has not yet been drawn ");
227 return;
228 }
230}
231
232////////////////////////////////////////////////////////////////////////////////
233/// Perform 2 times a gaussian fit to fetch the center of the histo.
234/// To get the second fit range get an interval that tries to keep into account
235/// the skewness of the distribution.
236
238 // Get the center of the histo
239
240 TString optfit = "Q0";
241 if (display_result) optfit = "Q";
242
243 TH1F* histo = static_cast<TH1F*>(histo_orig->Clone());
244
245 // get the histo x extremes
246 double x_min = histo->GetXaxis()->GetXmin();
247 double x_max = histo->GetXaxis()->GetXmax();
248
249 // First fit!
250
251 TF1* gauss = new TF1("mygauss", "gauss", x_min, x_max);
252
253 gauss->SetParameter("Constant",histo->GetEntries());
254 gauss->SetParameter("Mean",histo->GetMean());
255 gauss->SetParameter("Sigma",histo->GetRMS());
256
257 histo->Fit(gauss,optfit);
258
259 // Second fit!
260 double sigma = gauss->GetParameter("Sigma");
261 double mean = gauss->GetParameter("Mean");
262
263 delete gauss;
264
265 std::cout << "Center is 1st pass = " << mean << std::endl;
266
267 double skewness = histo->GetSkewness();
268
269 x_min = mean - n_rms*sigma - sigma*skewness/2;
270 x_max = mean + n_rms*sigma - sigma*skewness/2;
271
272 TF1* gauss2 = new TF1("mygauss2", "gauss", x_min, x_max);
273 gauss2->SetParameter("Mean",mean);
274
275 // second fit : likelihood fit
276 optfit += "L";
277 histo->Fit(gauss2,optfit,"", x_min, x_max);
278
279
280 double center = gauss2->GetParameter("Mean");
281
282 if (display_result) {
283 histo->Draw();
284 gauss2->Draw("same");
285 }
286 else {
287 delete histo;
288 }
289 delete gauss2;
290
291 return center;
292
293
294}
295
296////////////////////////////////////////////////////////////////////////////////
297/// We let an horizontal bar go down and we stop when we have the integral
298/// equal to the desired one.
299
300double* HybridPlot::GetHistoPvals (TH1* histo, double percentage){
301
302 if (percentage>1){
303 std::cerr << "Percentage greater or equal to 1!\n";
304 return nullptr;
305 }
306
307 // Get the integral of the histo
308 double* h_integral=histo->GetIntegral();
309
310 // Start the iteration
311 std::map<int,int> extremes_map;
312
313 for (int i=0;i<histo->GetNbinsX();++i){
314 for (int j=0;j<histo->GetNbinsX();++j){
315 double integral = h_integral[j]-h_integral[i];
316 if (integral>percentage){
317 extremes_map[i]=j;
318 break;
319 }
320 }
321 }
322
323 // Now select the couple of extremes which have the lower bin content diff
324 std::map<int,int>::iterator it;
325 int a;
326 int b;
327 double left_bin_center(0.);
328 double right_bin_center(0.);
329 double diff=10e40;
330 double current_diff;
331 for (it = extremes_map.begin();it != extremes_map.end();++it){
332 a=it->first;
333 b=it->second;
334 current_diff=std::abs(histo->GetBinContent(a)-histo->GetBinContent(b));
335 if (current_diff<diff){
336 //std::cout << "a=" << a << " b=" << b << std::endl;
340 }
341 }
342
343 double* d = new double[2];
346 return d;
347}
348
349////////////////////////////////////////////////////////////////////////////////
350/// Get the median of an histogram.
351
353
354 //int xbin_median;
355 double* integral = histo->GetIntegral();
356 int median_i = 0;
357 for (int j = 0; j < histo->GetNbinsX(); j++) {
358 if (integral[j]<0.5)
359 median_i = j;
360 }
361
362 double median_x =
363 histo->GetBinCenter(median_i)+
364 (histo->GetBinCenter(median_i+1)-
365 histo->GetBinCenter(median_i))*
366 (0.5-integral[median_i])/(integral[median_i+1]-integral[median_i]);
367
368 return median_x;
369}
bool fVerbose
The verbosity flag.
Definition HLFactory.h:70
#define d(i)
Definition RSha256.hxx:102
#define b(i)
Definition RSha256.hxx:100
#define a(i)
Definition RSha256.hxx:99
@ kRed
Definition Rtypes.h:66
@ kBlack
Definition Rtypes.h:65
@ kBlue
Definition Rtypes.h:66
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
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 filename
char name[80]
Definition TGX11.cxx:110
R__EXTERN TStyle * gStyle
Definition TStyle.h:442
#define gPad
const_iterator begin() const
const_iterator end() const
double GetMedian(TH1 *histo)
Get the median of an histogram.
TH1F * fB_histo_shaded
The b Histo shaded.
Definition HybridPlot.h:115
TH1F * fSb_histo_shaded
The sb Histo shaded.
Definition HybridPlot.h:113
TVirtualPad * fPad
The pad where it has been drawn.
Definition HybridPlot.h:118
HybridPlot(const char *name, const char *title, const std::vector< double > &sb_values, const std::vector< double > &b_values, double testStat_data, int n_bins, bool verbosity=true)
Constructor.
void Draw(const char *options="") override
Draw on current pad.
TH1F * fSb_histo
The sb Histo.
Definition HybridPlot.h:112
void DumpToImage(const char *filename)
Write an image on disk.
TLine * fData_testStat_line
The line for the data value of the test statistic.
Definition HybridPlot.h:116
TLegend * fLegend
The legend of the plot.
Definition HybridPlot.h:117
TH1F * fB_histo
The b Histo.
Definition HybridPlot.h:114
~HybridPlot() override
Destructor.
double * GetHistoPvals(TH1 *histo, double percentage)
Get the "effective sigmas" of the histo, call delete [] res to release memory.
double GetHistoCenter(TH1 *histo, double n_rms=1, bool display_result=false)
Get the center of the histo.
void DumpToFile(const char *RootFileName, const char *options)
All the objects are written to rootfile.
virtual void SetFillColor(Color_t fcolor)
Set the fill area color.
Definition TAttFill.h:38
virtual void SetFillStyle(Style_t fstyle)
Set the fill area style.
Definition TAttFill.h:40
virtual void SetLineWidth(Width_t lwidth)
Set the line width.
Definition TAttLine.h:45
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition TAttLine.h:42
Double_t GetXmax() const
Definition TAxis.h:142
Double_t GetXmin() const
Definition TAxis.h:141
1-Dim function class
Definition TF1.h:234
virtual void SetParameter(Int_t param, Double_t value)
Definition TF1.h:671
virtual Double_t GetParameter(Int_t ipar) const
Definition TF1.h:544
A ROOT file is an on-disk file, usually with extension .root, that stores objects in a file-system-li...
Definition TFile.h:131
1-D histogram with a float per channel (see TH1 documentation)
Definition TH1.h:647
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:59
virtual Double_t GetBinCenter(Int_t bin) const
Return bin center for 1D histogram.
Definition TH1.cxx:9133
void SetTitle(const char *title) override
Change/set the title.
Definition TH1.cxx:6721
virtual Double_t GetMean(Int_t axis=1) const
For axis = 1,2 or 3 returns the mean value of the histogram along X,Y or Z axis.
Definition TH1.cxx:7527
virtual Double_t GetSkewness(Int_t axis=1) const
Definition TH1.cxx:7663
virtual TH1 * DrawNormalized(Option_t *option="", Double_t norm=1) const
Draw a normalized copy of this histogram.
Definition TH1.cxx:3116
TAxis * GetXaxis()
Definition TH1.h:341
virtual TFitResultPtr Fit(const char *formula, Option_t *option="", Option_t *goption="", Double_t xmin=0, Double_t xmax=0)
Fit histogram with function fname.
Definition TH1.cxx:3876
virtual Double_t GetMaximum(Double_t maxval=FLT_MAX) const
Return maximum value smaller than maxval of bins in the range, unless the value has been overridden b...
Definition TH1.cxx:8537
virtual Int_t GetNbinsX() const
Definition TH1.h:314
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition TH1.cxx:3316
void Draw(Option_t *option="") override
Draw this histogram with options.
Definition TH1.cxx:3038
Double_t GetRMS(Int_t axis=1) const
This function returns the Standard Deviation (Sigma) of the distribution not the Root Mean Square (RM...
Definition TH1.h:336
virtual void SetBinContent(Int_t bin, Double_t content)
Set bin content see convention for numbering bins in TH1::GetBin In case the bin number is greater th...
Definition TH1.cxx:9214
virtual Double_t GetEntries() const
Return the current number of entries.
Definition TH1.cxx:4401
virtual Double_t * GetIntegral()
Return a pointer to the array of bins integral.
Definition TH1.cxx:2558
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition TH1.cxx:5064
TObject * Clone(const char *newname="") const override
Make a complete copy of the underlying object.
Definition TH1.cxx:2724
virtual Double_t GetSumOfWeights() const
Return the sum of weights excluding under/overflows.
Definition TH1.cxx:7909
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:320
void Draw(Option_t *option="") override
Draw this legend with its current attributes.
Definition TLegend.cxx:425
Use the TLine constructor to create a simple line.
Definition TLine.h:22
Double_t GetX1() const
Definition TLine.h:50
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:174
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
virtual Int_t Write(const char *name=nullptr, Int_t option=0, Int_t bufsize=0)
Write this object to the current directory.
Definition TObject.cxx:947
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:1054
virtual void Draw(Option_t *option="")
Default Draw method for all objects.
Definition TObject.cxx:293
virtual void SetName(const char *name="")
Definition TPave.h:81
Basic string class.
Definition TString.h:139
void SetOptStat(Int_t stat=1)
The type of information printed in the histogram statistics box can be selected via the parameter mod...
Definition TStyle.cxx:1642
void Print(const char *filename="") const override=0
This method must be overridden when a class wants to print itself.
const Double_t sigma
Namespace for the RooStats classes.
Definition CodegenImpl.h:58