00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016 #include "assert.h"
00017 #include <cmath>
00018 #include <iostream>
00019 #include <map>
00020
00021 #include "RooStats/HybridPlot.h"
00022 #include "TStyle.h"
00023 #include "TF1.h"
00024 #include "TAxis.h"
00025 #include "TH1.h"
00026 #include "TLine.h"
00027 #include "TLegend.h"
00028 #include "TFile.h"
00029 #include "TVirtualPad.h"
00030
00031 #include <algorithm>
00032
00033
00034 ClassImp(RooStats::HybridPlot)
00035
00036 using namespace RooStats;
00037
00038
00039
00040 HybridPlot::HybridPlot(const char* name,
00041 const char* title,
00042 const std::vector<double> & sb_vals,
00043 const std::vector<double> & b_vals,
00044 double testStat_data,
00045 int n_bins,
00046 bool verbosity):
00047 TNamed(name,title),
00048 fSb_histo(NULL),
00049 fSb_histo_shaded(NULL),
00050 fB_histo(NULL),
00051 fB_histo_shaded(NULL),
00052 fData_testStat_line(0),
00053 fLegend(0),
00054 fPad(0),
00055 fVerbose(verbosity)
00056 {
00057
00058
00059 int nToysSB = sb_vals.size();
00060 int nToysB = sb_vals.size();
00061 assert (nToysSB >0);
00062 assert (nToysB >0);
00063
00064
00065 double min = *std::min_element(sb_vals.begin(), sb_vals.end());
00066 double max = *std::max_element(sb_vals.begin(), sb_vals.end());
00067
00068 double min_b = *std::min_element(b_vals.begin(), b_vals.end());
00069 double max_b = *std::max_element(b_vals.begin(), b_vals.end());
00070
00071
00072 if ( min_b < min) min = min_b;
00073 if ( max_b > max) max = max_b;
00074
00075 if (testStat_data<min) min = testStat_data;
00076 if (testStat_data>max) max = testStat_data;
00077
00078 min *= 1.1;
00079 max *= 1.1;
00080
00081
00082
00083 fSb_histo = new TH1F ("SB_model",title,n_bins,min,max);
00084 fSb_histo->SetTitle(fSb_histo->GetTitle());
00085 fSb_histo->SetLineColor(kBlue);
00086 fSb_histo->GetXaxis()->SetTitle("test statistics");
00087 fSb_histo->SetLineWidth(2);
00088
00089 fB_histo = new TH1F ("B_model",title,n_bins,min,max);
00090 fB_histo->SetTitle(fB_histo->GetTitle());
00091 fB_histo->SetLineColor(kRed);
00092 fB_histo->GetXaxis()->SetTitle("test statistics");
00093 fB_histo->SetLineWidth(2);
00094
00095 for (int i=0;i<nToysSB;++i) fSb_histo->Fill(sb_vals[i]);
00096 for (int i=0;i<nToysB;++i) fB_histo->Fill(b_vals[i]);
00097
00098 double histos_max_y = fSb_histo->GetMaximum();
00099 double line_hight = histos_max_y/nToysSB;
00100 if (histos_max_y<fB_histo->GetMaximum()) histos_max_y = fB_histo->GetMaximum()/nToysB;
00101
00102
00103 fData_testStat_line = new TLine(testStat_data,0,testStat_data,line_hight);
00104 fData_testStat_line->SetLineWidth(3);
00105 fData_testStat_line->SetLineColor(kBlack);
00106
00107
00108 double golden_section = (std::sqrt(5.)-1)/2;
00109
00110 fLegend = new TLegend(0.75,0.95-0.2*golden_section,0.95,0.95);
00111 TString title_leg="test statistics distributions ";
00112 title_leg+=sb_vals.size();
00113 title_leg+=" toys";
00114 fLegend->SetName(title_leg.Data());
00115 fLegend->AddEntry(fSb_histo,"SB toy datasets");
00116 fLegend->AddEntry(fB_histo,"B toy datasets");
00117 fLegend->AddEntry((TLine*)fData_testStat_line,"test statistics on data","L");
00118 fLegend->SetFillColor(0);
00119
00120 }
00121
00122
00123
00124 HybridPlot::~HybridPlot()
00125 {
00126
00127
00128 if (fSb_histo) delete fSb_histo;
00129 if (fB_histo) delete fB_histo;
00130 if (fSb_histo_shaded) delete fSb_histo_shaded;
00131 if (fB_histo_shaded) delete fB_histo_shaded;
00132 if (fData_testStat_line) delete fData_testStat_line;
00133 if (fLegend) delete fLegend;
00134 }
00135
00136
00137
00138 void HybridPlot::Draw(const char* )
00139 {
00140
00141
00142
00143
00144 gStyle->SetOptStat(0);
00145
00146
00147
00148 if (fSb_histo->GetMaximum()>fB_histo->GetMaximum()){
00149 fSb_histo->DrawNormalized();
00150 fB_histo->DrawNormalized("same");
00151 }
00152 else{
00153 fB_histo->DrawNormalized();
00154 fSb_histo->DrawNormalized("same");
00155 }
00156
00157
00158 fB_histo_shaded = (TH1F*)fB_histo->Clone("b_shaded");
00159 fB_histo_shaded->SetFillStyle(3005);
00160 fB_histo_shaded->SetFillColor(kRed);
00161
00162 fSb_histo_shaded = (TH1F*)fSb_histo->Clone("sb_shaded");
00163 fSb_histo_shaded->SetFillStyle(3004);
00164 fSb_histo_shaded->SetFillColor(kBlue);
00165
00166
00167 double data_m2lnq= fData_testStat_line->GetX1();
00168 for (int i=1;i<=fSb_histo->GetNbinsX();++i){
00169 if (fSb_histo->GetBinCenter(i)<data_m2lnq){
00170 fSb_histo_shaded->SetBinContent(i,0);
00171 fB_histo_shaded->SetBinContent(i,fB_histo->GetBinContent(i)/fB_histo->GetSumOfWeights());
00172 }
00173 else{
00174 fB_histo_shaded->SetBinContent(i,0);
00175 fSb_histo_shaded->SetBinContent(i,fSb_histo->GetBinContent(i)/fSb_histo->GetSumOfWeights());
00176 }
00177 }
00178
00179
00180 fB_histo_shaded->Draw("same");
00181 fSb_histo_shaded->Draw("same");
00182
00183
00184 fData_testStat_line->Draw("same");
00185
00186
00187 fLegend->Draw("same");
00188
00189 if (gPad) {
00190 gPad->SetName(GetName());
00191 gPad->SetTitle(GetTitle() );
00192 }
00193
00194 fPad = gPad;
00195
00196 }
00197
00198
00199
00200 void HybridPlot::DumpToFile (const char* RootFileName, const char* options)
00201 {
00202
00203
00204 TFile ofile(RootFileName,options);
00205 ofile.cd();
00206
00207
00208 fSb_histo->Write();
00209 fB_histo->Write();
00210
00211
00212 if (fB_histo_shaded!=NULL && fSb_histo_shaded!=NULL){
00213 fB_histo_shaded->Write();
00214 fSb_histo_shaded->Write();
00215 }
00216
00217
00218 fData_testStat_line->Write("Measured test statistics line tag");
00219
00220
00221 fLegend->Write();
00222
00223 ofile.Close();
00224
00225 }
00226
00227 void HybridPlot::DumpToImage(const char * filename) {
00228 if (!fPad) {
00229 Error("HybridPlot","Hybrid plot has not yet been drawn ");
00230 return;
00231 }
00232 fPad->Print(filename);
00233 }
00234
00235
00236
00237
00238
00239
00240
00241
00242
00243
00244
00245 double HybridPlot::GetHistoCenter(TH1* histo_orig, double n_rms, bool display_result){
00246
00247
00248 TString optfit = "Q0";
00249 if (display_result) optfit = "Q";
00250
00251 TH1F* histo = (TH1F*)histo_orig->Clone();
00252
00253
00254 double x_min = histo->GetXaxis()->GetXmin();
00255 double x_max = histo->GetXaxis()->GetXmax();
00256
00257
00258
00259 TF1* gaus = new TF1("mygaus", "gaus", x_min, x_max);
00260
00261 gaus->SetParameter("Constant",histo->GetEntries());
00262 gaus->SetParameter("Mean",histo->GetMean());
00263 gaus->SetParameter("Sigma",histo->GetRMS());
00264
00265 histo->Fit(gaus,optfit);
00266
00267
00268 double sigma = gaus->GetParameter("Sigma");
00269 double mean = gaus->GetParameter("Mean");
00270
00271 delete gaus;
00272
00273 std::cout << "Center is 1st pass = " << mean << std::endl;
00274
00275 double skewness = histo->GetSkewness();
00276
00277 x_min = mean - n_rms*sigma - sigma*skewness/2;
00278 x_max = mean + n_rms*sigma - sigma*skewness/2;;
00279
00280 TF1* gaus2 = new TF1("mygaus2", "gaus", x_min, x_max);
00281 gaus2->SetParameter("Mean",mean);
00282
00283
00284 optfit += "L";
00285 histo->Fit(gaus2,optfit,"", x_min, x_max);
00286
00287
00288 double center = gaus2->GetParameter("Mean");
00289
00290 if (display_result) {
00291 histo->Draw();
00292 gaus2->Draw("same");
00293 }
00294 else {
00295 delete histo;
00296 }
00297 delete gaus2;
00298
00299 return center;
00300
00301
00302 }
00303
00304
00305
00306
00307
00308
00309 double* HybridPlot::GetHistoPvals (TH1* histo, double percentage){
00310
00311 if (percentage>1){
00312 std::cerr << "Percentage greater or equal to 1!\n";
00313 return NULL;
00314 }
00315
00316
00317 double* h_integral=histo->GetIntegral();
00318
00319
00320 std::map<int,int> extremes_map;
00321
00322 for (int i=0;i<histo->GetNbinsX();++i){
00323 for (int j=0;j<histo->GetNbinsX();++j){
00324 double integral = h_integral[j]-h_integral[i];
00325 if (integral>percentage){
00326 extremes_map[i]=j;
00327 break;
00328 }
00329 }
00330 }
00331
00332
00333 std::map<int,int>::iterator it;
00334 int a,b;
00335 double left_bin_center(0.),right_bin_center(0.);
00336 double diff=10e40;
00337 double current_diff;
00338 for (it = extremes_map.begin();it != extremes_map.end();++it){
00339 a=it->first;
00340 b=it->second;
00341 current_diff=std::fabs(histo->GetBinContent(a)-histo->GetBinContent(b));
00342 if (current_diff<diff){
00343
00344 diff=current_diff;
00345 left_bin_center=histo->GetBinCenter(a);
00346 right_bin_center=histo->GetBinCenter(b);
00347 }
00348 }
00349
00350 double* d = new double[2];
00351 d[0]=left_bin_center;
00352 d[1]=right_bin_center;
00353 return d;
00354 }
00355
00356
00357
00358
00359
00360 double HybridPlot::GetMedian(TH1* histo){
00361
00362
00363 Double_t* integral = histo->GetIntegral();
00364 int median_i = 0;
00365 for (int j=0;j<histo->GetNbinsX(); j++)
00366 if (integral[j]<0.5)
00367 median_i = j;
00368
00369 double median_x =
00370 histo->GetBinCenter(median_i)+
00371 (histo->GetBinCenter(median_i+1)-
00372 histo->GetBinCenter(median_i))*
00373 (0.5-integral[median_i])/(integral[median_i+1]-integral[median_i]);
00374
00375 return median_x;
00376 }
00377
00378