CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_4_5_patch3/src/TopQuarkAnalysis/TopTools/src/LRHelpFunctions.cc

Go to the documentation of this file.
00001 #include "TopQuarkAnalysis/TopTools/interface/LRHelpFunctions.h"
00002 #include "TopQuarkAnalysis/TopTools/test/tdrstyle.C"
00003 
00004 // constructors
00005 LRHelpFunctions::LRHelpFunctions() {
00006   constructPurity = false;
00007   setTDRStyle();
00008   gStyle->SetCanvasDefW(900);
00009 }
00010 
00011 
00012 LRHelpFunctions::LRHelpFunctions(std::vector<int> obsNr, int nrBins, std::vector<double> obsMin, std::vector<double> obsMax, std::vector<const char*> functions) { 
00013   obsNumbers = obsNr;
00014   constructPurity = false;
00015   setTDRStyle();
00016   gStyle->SetCanvasDefW(900);
00017   for(size_t o=0; o<obsNr.size(); o++){
00018     // create Signal, background and s/(s+b) histogram
00019     TString htS  = "Obs"; htS  += obsNr[o]; htS += "_S"; 
00020     TString htB  = "Obs"; htB  += obsNr[o]; htB += "_B"; 
00021     TString htSB = "Obs"; htSB += obsNr[o]; htSB += "_SoverSplusB"; 
00022     hObsS.push_back( new TH1F(htS,htS,nrBins,obsMin[o],obsMax[o]) );
00023     hObsB.push_back( new TH1F(htB,htB,nrBins,obsMin[o],obsMax[o]) );
00024     hObsSoverSplusB.push_back( new TH1F(htSB,htSB,nrBins,obsMin[o],obsMax[o]) );
00025     
00026     // create the correlation 2D plots among the observables (for signal)
00027     for (size_t o2 = o+1; o2 < obsNr.size(); o2++) {
00028       TString hcorr  = "Corr_Obs"; hcorr  += obsNr[o]; hcorr += "_Obs";  hcorr  += obsNr[o2]; 
00029       hObsCorr.push_back( new TH2F(hcorr,hcorr,nrBins,obsMin[o],obsMax[o],nrBins,obsMin[o2],obsMax[o2]) );
00030     }    
00031     
00032     // create fit functions
00033     TString ftSB = "F_Obs"; ftSB += obsNr[o]; ftSB += "_SoverSplusB"; 
00034     fObsSoverSplusB.push_back( new TF1(ftSB,functions[o],hObsS[o]->GetXaxis()->GetXmin(),hObsS[o]->GetXaxis()->GetXmax()) );
00035   }
00036 }
00037 
00038 void LRHelpFunctions::recreateFitFct(std::vector<int> obsNr, std::vector<const char*> functions) {
00039   if (!fObsSoverSplusB.empty()) fObsSoverSplusB.clear();
00040   for(size_t o=0; o<obsNr.size(); o++){
00041     // create fit functions
00042     TString ftSB = "F_Obs"; ftSB += obsNr[o]; ftSB += "_SoverSplusB";
00043     fObsSoverSplusB.push_back( new TF1(ftSB,functions[o],hObsS[o]->GetXaxis()->GetXmin(),hObsS[o]->GetXaxis()->GetXmax()) );
00044   }
00045 }
00046 
00047 LRHelpFunctions::LRHelpFunctions(int nrLRbins, double LRmin, double LRmax, const char* LRfunction) { 
00048   constructPurity = true;
00049   setTDRStyle();
00050   gStyle->SetCanvasDefW(900);
00051 
00052   // create LR histograms
00053   hLRtotS = new TH1F("hLRtotS","hLRtotS",nrLRbins,LRmin,LRmax);
00054   hLRtotS->GetXaxis()->SetTitle("Combined LR");
00055   hLRtotB = new TH1F("hLRtotB","hLRtotB",nrLRbins,LRmin,LRmax);
00056   hLRtotB->GetXaxis()->SetTitle("Combined LR");
00057   hLRtotSoverSplusB = new TH1F("hLRtotSoverSplusB","hLRtotSoverSplusB",nrLRbins,LRmin,LRmax);
00058     
00059   // create LR fit function
00060   fLRtotSoverSplusB = new TF1("fLRtotSoverSplusB",LRfunction,LRmin,LRmax);
00061 }
00062 
00063 
00064 // destructor
00065 LRHelpFunctions::~LRHelpFunctions() {}
00066 
00067 
00068 // member function to initialize the LR hists and fits
00069 void LRHelpFunctions::initLRHistsAndFits(int nrLRbins, double LRmin, double LRmax, const char* LRfunction){ 
00070   constructPurity = true;
00071   // create LR histograms
00072   hLRtotS = new TH1F("hLRtotS","hLRtotS",nrLRbins,LRmin,LRmax);
00073   hLRtotB = new TH1F("hLRtotB","hLRtotB",nrLRbins,LRmin,LRmax);
00074   hLRtotSoverSplusB = new TH1F("hLRtotSoverSplusB","hLRtotSoverSplusB",nrLRbins,LRmin,LRmax);    
00075   // create LR fit function
00076   fLRtotSoverSplusB = new TF1("fLRtotSoverSplusB",LRfunction,LRmin,LRmax);
00077 }
00078 
00079 
00080 
00081 // member function to set initial values to the observable fit function
00082 void LRHelpFunctions::setObsFitParameters(int obs,std::vector<double> fitPars){
00083   for(size_t fit=0; fit<fObsSoverSplusB.size(); fit++){
00084     TString fn = "_Obs"; fn += obs;
00085     if(((TString)fObsSoverSplusB[fit]->GetName()).Contains(fn)){
00086       for(size_t p=0; p<fitPars.size(); p++){
00087         fObsSoverSplusB[fit]->SetParameter(p,fitPars[p]);
00088       }
00089     }
00090   }
00091 }
00092 
00093 
00094 
00095 
00096 
00097 // member function to add observable values to the signal histograms
00098 void LRHelpFunctions::fillToSignalHists(std::vector<double> obsVals, double weight){
00099   int hIndex = 0;
00100   for(size_t o=0; o<obsVals.size(); o++) {
00101     hObsS[o]->Fill(obsVals[o], weight);
00102     for(size_t o2=o+1; o2<obsVals.size(); o2++) {
00103       hObsCorr[hIndex] -> Fill(obsVals[o],obsVals[o2], weight);
00104       ++hIndex;
00105     }
00106   }
00107 }
00108 
00109 // member function to add observable values to the signal histograms
00110 void LRHelpFunctions::fillToSignalHists(int obsNbr, double obsVal, double weight){
00111   TString obs = "Obs"; obs += obsNbr; obs += "_";
00112   for(size_t f = 0; f<hObsS.size(); f++){
00113     if(((TString)(hObsS[f]->GetName())).Contains(obs)) {
00114       hObsS[f]->Fill(obsVal, weight);
00115       return;
00116     }
00117   }
00118 }
00119 
00120 // member function to add observable values to the signal histograms
00121 void LRHelpFunctions::fillToSignalCorrelation(int obsNbr1, double obsVal1, int obsNbr2,
00122         double obsVal2, double weight){
00123   TString hcorr  = "Corr_Obs"; hcorr  += std::min(obsNbr1,obsNbr2) ; hcorr += "_Obs";  hcorr  += std::max(obsNbr1, obsNbr2);
00124   for(size_t i=0; i<hObsCorr.size(); i++) {
00125     if(((TString)(hObsCorr[i]->GetName())) == hcorr) {
00126       if (obsNbr1 < obsNbr2) {
00127         hObsCorr[i] -> Fill(obsVal1,obsVal2, weight);
00128       } else {
00129         hObsCorr[i] -> Fill(obsVal2,obsVal1, weight);
00130       }
00131       return;
00132     }
00133   }
00134 }
00135 
00136 
00137 // member function to add observable values to the background histograms
00138 void LRHelpFunctions::fillToBackgroundHists(std::vector<double> obsVals, double weight){
00139   for(size_t o=0; o<obsVals.size(); o++) hObsB[o]->Fill(obsVals[o], weight);
00140 }
00141 
00142 // member function to add observable values to the signal histograms
00143 void LRHelpFunctions::fillToBackgroundHists(int obsNbr, double obsVal, double weight){
00144   TString obs = "Obs"; obs += obsNbr; obs += "_";
00145   for(size_t f = 0; f<hObsB.size(); f++){
00146     if(((TString)(hObsB[f]->GetName())).Contains(obs)) {
00147       hObsB[f]->Fill(obsVal, weight);
00148       return;
00149     }
00150   }
00151 }
00152 
00153 // member function to normalize the S and B spectra
00154 void LRHelpFunctions::normalizeSandBhists(){
00155   for(size_t o=0; o<hObsS.size(); o++){
00156     // count entries in each histo. Do it this way instead of GetEntries method,
00157     // since the latter does not account for the weights.
00158     double nrSignEntries = 0, nrBackEntries = 0;
00159     for (int i=0; i<=hObsS[o]->GetNbinsX()+1; i++){
00160       nrSignEntries += hObsS[o]->GetBinContent(i);
00161       nrBackEntries += hObsB[o]->GetBinContent(i);
00162      }
00163     for(int b=0; b <= hObsS[o]->GetNbinsX()+1; b++){
00164       hObsS[o]->SetBinContent(b,hObsS[o]->GetBinContent(b)/(nrSignEntries));
00165       hObsB[o]->SetBinContent(b,hObsB[o]->GetBinContent(b)/(nrBackEntries));
00166       hObsS[o]->SetBinError(b,hObsS[o]->GetBinError(b)/(nrSignEntries));
00167       hObsB[o]->SetBinError(b,hObsB[o]->GetBinError(b)/(nrBackEntries));
00168     }
00169     //std::cout<<"Integral for obs"<<o<<" S: "<<hObsS[o]->Integral(0,10000)<<" & obs"<<o<<" B: "<<hObsB[o]->Integral(0,10000)<<std::endl;
00170   }
00171 }
00172 
00173  
00174     
00175 // member function to produce and fit the S/(S+B) histograms
00176 void LRHelpFunctions::makeAndFitSoverSplusBHists(){
00177   for(size_t o=0; o<hObsS.size(); o++){
00178     for(int b=0; b <= hObsS[o]->GetNbinsX()+1; b++){
00179       if ((hObsS[o]->GetBinContent(b)+ hObsB[o]->GetBinContent(b)) > 0) {
00180         hObsSoverSplusB[o]->SetBinContent(b, hObsS[o]->GetBinContent(b) / (hObsS[o]->GetBinContent(b) + hObsB[o]->GetBinContent(b)));
00181         double error = sqrt( pow(hObsS[o]->GetBinError(b) * hObsB[o]->GetBinContent(b)/pow(hObsS[o]->GetBinContent(b) + hObsB[o]->GetBinContent(b),2),2)
00182                            + pow(hObsB[o]->GetBinError(b) * hObsS[o]->GetBinContent(b)/pow(hObsS[o]->GetBinContent(b) + hObsB[o]->GetBinContent(b),2),2) );
00183         hObsSoverSplusB[o]->SetBinError(b,error);
00184       }
00185     }
00186     hObsSoverSplusB[o]->Fit(fObsSoverSplusB[o]->GetName(),"RQ");
00187   }
00188 }
00189 
00190 
00191 // member function to read the observable hists & fits from a root-file
00192 void LRHelpFunctions::readObsHistsAndFits(TString fileName, std::vector<int> observables, bool readLRplots){
00193   hObsS.clear();
00194   hObsB.clear();
00195   hObsSoverSplusB.clear();
00196   fObsSoverSplusB.clear();
00197   TFile *fitFile = new TFile(fileName, "READ");
00198   if(observables[0] == -1){  
00199     std::cout<<" ... will read hists and fit for all available observables in file "<<fileName<<std::endl;
00200     TList *list = fitFile->GetListOfKeys();
00201     TIter next(list);
00202     TKey *el;
00203     while ((el = (TKey*)next())) {
00204       TString keyName = el->GetName();
00205       if(keyName.Contains("F_") && keyName.Contains("_SoverSplusB")){
00206         fObsSoverSplusB.push_back( new TF1(*((TF1*) el -> ReadObj())));
00207       }
00208       else if(keyName.Contains("_SoverSplusB")){
00209         hObsSoverSplusB.push_back( new TH1F(*((TH1F*) el -> ReadObj())));
00210       }
00211       else if(keyName.Contains("_S")){
00212         hObsS.push_back( new TH1F(*((TH1F*) el -> ReadObj())));
00213       }
00214       else if(keyName.Contains("_B")){
00215         hObsB.push_back( new TH1F(*((TH1F*) el -> ReadObj())));
00216       }
00217       else if(keyName.Contains("Corr")){
00218         hObsCorr.push_back( new TH2F(*((TH2F*) el -> ReadObj())));
00219       }
00220     }
00221   }
00222   else{
00223     obsNumbers = observables;
00224     for(unsigned int obs = 0; obs < observables.size(); obs++){
00225       std::cout<<"  ... will read hists and fit for obs "<<observables[obs]<<" from file "<<fileName<<std::endl;
00226       TString hStitle  = "Obs";   hStitle  += observables[obs];   hStitle += "_S";
00227       hObsS.push_back( new TH1F(*((TH1F*)fitFile->GetKey(hStitle)->ReadObj())));
00228       TString hBtitle  = "Obs";   hBtitle  += observables[obs];   hBtitle += "_B";
00229       hObsB.push_back( new TH1F(*((TH1F*)fitFile->GetKey(hBtitle)->ReadObj())));
00230       TString hSBtitle = "Obs";   hSBtitle += observables[obs];   hSBtitle += "_SoverSplusB";
00231       TString fSBtitle = "F_";  fSBtitle += hSBtitle;
00232       hObsSoverSplusB.push_back( new TH1F(*((TH1F*)fitFile->GetKey(hSBtitle)->ReadObj())));
00233       fObsSoverSplusB.push_back( new TF1(*((TF1*)fitFile->GetKey(fSBtitle)->ReadObj())));
00234       for(unsigned int obs2 = obs+1; obs2 < observables.size(); obs2++){
00235         TString hCorrtitle  = "Corr_Obs"; hCorrtitle  += observables[obs];   
00236                 hCorrtitle += "_Obs";     hCorrtitle  += observables[obs2]; 
00237         hObsCorr.push_back( new TH2F(*((TH2F*)fitFile->GetKey(hCorrtitle)->ReadObj())));
00238       }
00239     }
00240   }
00241   
00242   if(readLRplots){
00243     constructPurity = true;
00244     std::cout<<"  ... will LR s and B histos from file "<<fileName<<std::endl;
00245     hLRtotS = new TH1F(*((TH1F*)fitFile->GetKey("hLRtotS")->ReadObj()));
00246     hLRtotB = new TH1F(*((TH1F*)fitFile->GetKey("hLRtotB")->ReadObj()));
00247     hLRtotSoverSplusB = new TH1F(*((TH1F*)fitFile->GetKey("hLRtotSoverSplusB")->ReadObj()));
00248     fLRtotSoverSplusB = new TF1(*((TF1*)fitFile->GetKey("fLRtotSoverSplusB") -> ReadObj()));
00249   }
00250 }
00251 
00252 
00253 
00254 
00255 // member function to store all observable plots and fits to a ROOT-file
00256 void  LRHelpFunctions::storeToROOTfile(TString fname){
00257   TFile fOut(fname,"RECREATE");
00258   fOut.cd();
00259   for(size_t o=0; o<hObsS.size(); o++){
00260     hObsS[o]            -> Write();
00261     hObsB[o]            -> Write();
00262     hObsSoverSplusB[o]  -> Write();
00263     fObsSoverSplusB[o]  -> Write();
00264   }
00265   int hIndex = 0;
00266   for(size_t o=0; o<hObsS.size(); o++) {
00267     for(size_t o2=o+1; o2<hObsS.size(); o2++) {
00268       hObsCorr[hIndex] -> Write();
00269       ++ hIndex;
00270     }
00271   }
00272   if(constructPurity){
00273     hLRtotS             -> Write();
00274     hLRtotB             -> Write();
00275     hLRtotSoverSplusB   -> Write();
00276     fLRtotSoverSplusB   -> Write();
00277     hEffvsPur           -> Write();
00278     hLRValvsPur         -> Write();
00279     hLRValvsEff         -> Write();
00280   }
00281   fOut.cd();
00282   fOut.Write();
00283   fOut.Close();
00284 }
00285 
00286 
00287 
00288 
00289 
00290 // member function to make some simple control plots and store them in a ps-file
00291 void  LRHelpFunctions::storeControlPlots(TString fname){  
00292   TCanvas c("dummy","",1);
00293   c.Print(fname + "[","landscape");
00294   for(size_t o=0; o<hObsS.size(); o++) {
00295      TCanvas c2("c2","",1);
00296      c2.Divide(2,1);
00297      c2.cd(1);
00298      hObsS[o] -> SetLineColor(2);
00299      if(hObsB[o]->GetMaximum() > hObsS[o]->GetMaximum()){
00300        hObsB[o] -> Draw("hist");
00301        hObsS[o] -> Draw("histsame");
00302      }
00303      else
00304      {
00305        hObsS[o] -> Draw("hist");
00306        hObsB[o] -> Draw("histsame");
00307      }
00308      c2.cd(2);
00309      hObsSoverSplusB[o] -> Draw();
00310      fObsSoverSplusB[o] -> Draw("same");
00311      c2.Print(fname,"Landscape");
00312   }
00313   
00314   int hIndex = 0;
00315   for(size_t o=0; o<hObsS.size(); o++) {
00316     for(size_t o2=o+1; o2<hObsS.size(); o2++) {
00317       TCanvas cc("cc","",1);
00318       hObsCorr[hIndex] -> Draw("box");
00319       TPaveText pt(0.5,0.87,0.98,0.93,"blNDC");
00320       pt.SetFillStyle(1);
00321       pt.SetFillColor(0);
00322       pt.SetBorderSize(0);
00323       TString tcorr = "Corr. of "; tcorr += (int)(100.*hObsCorr[hIndex] -> GetCorrelationFactor()); tcorr += " %";
00324       //TText *text = pt.AddText(tcorr);
00325       pt.Draw("same");
00326       ++hIndex;
00327       cc.Print(fname,"Landscape");
00328     }
00329   }
00330   
00331   if(constructPurity){
00332     TCanvas c3("c3","",1);
00333     c3.Divide(2,1);
00334     c3.cd(1);
00335     hLRtotB -> Draw();
00336     hLRtotS -> SetLineColor(2);
00337     hLRtotS -> Draw("same");
00338     c3.cd(2);
00339     hLRtotSoverSplusB -> Draw();
00340     c3.Print(fname,"Landscape");
00341   
00342     TCanvas c4("c4","",1);
00343     hEffvsPur -> Draw("AL*");
00344     c4.Print(fname,"Landscape");
00345 
00346     TCanvas c5("c5","",1);
00347     hLRValvsPur -> Draw("AL*");
00348     c5.Print(fname,"Landscape");
00349 
00350     TCanvas c6("c6","",1);
00351     hLRValvsEff -> Draw("AL*");
00352     c6.Print(fname,"Landscape");
00353   } 
00354   c.Print(fname + "]","landscape");
00355 }
00356    
00357 
00358 
00359 
00360 // member function to calculate the LR value, using the S/N definition
00361 double  LRHelpFunctions::calcLRval(std::vector<double> vals){
00362   double logLR = 0.;
00363   for(size_t o=0; o<fObsSoverSplusB.size(); o++){
00364     double SoverSplusN = fObsSoverSplusB[o]->Eval(vals[o]);
00365     double SoverN = 1./((1./SoverSplusN) - 1.);
00366     logLR += log(SoverN);
00367   }  
00368   return logLR;
00369 }
00370 
00371 
00372 // member function to calculate the LR value, using the definition that was
00373 // used in the P-TDR: S/(S+N)
00374  
00375 double LRHelpFunctions::calcPtdrLRval(std::vector<double> vals, bool useCorrelation) {
00376   double logLR = 1.;
00377   for(size_t o=0; o<fObsSoverSplusB.size(); o++){
00378     double SoverSplusN = fObsSoverSplusB[o]->Eval(vals[o]);
00379     if (SoverSplusN<0.0001) SoverSplusN=0.0001;
00380     if (useCorrelation){
00381       double corr = 0;
00382       for(size_t j=0; j<fObsSoverSplusB.size(); j++){
00383         if (o==j) ++corr;
00384         else {
00385           TString hcorr  = "Corr_Obs"; hcorr  += std::min(obsNumbers[o],obsNumbers[j]) ; hcorr += "_Obs";  hcorr  += std::max(obsNumbers[o],obsNumbers[j]);
00386           for(size_t i=0; i<hObsCorr.size(); i++) {
00387             if(((TString)(hObsCorr[i]->GetName())) == hcorr)
00388               corr += fabs(hObsCorr[i]->GetCorrelationFactor());
00389           }
00390         }
00391       }
00392      logLR *= pow(SoverSplusN/fObsSoverSplusB[o]->GetMaximum(), 1./corr);
00393     }
00394     else logLR *= SoverSplusN/fObsSoverSplusB[o]->GetMaximum();
00395   }
00396   //std::cout <<logLR<<std::endl;
00397   return logLR;
00398 }
00399 
00400 
00401 // member function to fill a signal contribution to the LR histogram
00402 void  LRHelpFunctions::fillLRSignalHist(double val, double weight){
00403   //  std::cout<< "@@@===> LRHelpf Signal LRval = "<< val<<std::endl;
00404   hLRtotS -> Fill(val, weight);
00405 }
00406 
00407 
00408 
00409 // member function to fill a background contribution to the LR histogram
00410 void  LRHelpFunctions::fillLRBackgroundHist(double val, double weight){
00411   // std::cout<< "@@@===> LRHelpf Backgroud LRval = "<< val<<std::endl;
00412   hLRtotB -> Fill(val, weight);
00413 }
00414   
00415   
00416    
00417 // member function to make and fit the purity vs. LRval histogram, and the purity vs. efficiency graph
00418 void  LRHelpFunctions::makeAndFitPurityHists(){
00419 
00420   for(int b=0; b<=hLRtotS->GetNbinsX(); b++) {
00421     float Sint = hLRtotS->Integral(b, hLRtotS->GetNbinsX()+1);
00422     float Bint = hLRtotB->Integral(b, hLRtotB->GetNbinsX()+1);
00423     if (Sint + Bint > 0) {
00424       hLRtotSoverSplusB->SetBinContent(b,1. * Sint / (Sint + Bint));
00425       hLRtotSoverSplusB->SetBinError(b,sqrt((pow(Sint*sqrt(Bint),2)+pow(Bint*sqrt(Sint),2)))/pow((Sint+Bint),2));
00426     }
00427   }
00428 
00429   hLRtotS->GetXaxis()->SetTitle("Combined LR ratio");
00430   hLRtotB->GetXaxis()->SetTitle("Combined LR ratio");
00431   hLRtotSoverSplusB->GetXaxis()->SetTitle("Cut on Combined LR");
00432   hLRtotSoverSplusB->GetYaxis()->SetTitle("Purity");
00433 
00434   hLRtotSoverSplusB -> Fit(fLRtotSoverSplusB->GetName(),"RQ");  
00435   double totSignal = hLRtotS->Integral(0,hLRtotS->GetNbinsX()+1);
00436   double Eff[200], Pur[200], LRVal[200];
00437   if (hLRtotS->GetNbinsX()>200) {
00438     std::cout << "Number of bins of LR histograms can not execeed 200!"<<std::endl;
00439     return;
00440   }
00441   for(int cut=0; (cut<=hLRtotS->GetNbinsX())&&(cut<200) ; cut ++){
00442         double LRcutVal = hLRtotS->GetBinLowEdge(cut);
00443         Eff[cut]   = hLRtotS->Integral(cut,hLRtotS->GetNbinsX()+1)/totSignal;
00444         Pur[cut]   = fLRtotSoverSplusB->Eval(LRcutVal);
00445         LRVal[cut] = LRcutVal;
00446   }
00447   hEffvsPur = new TGraph(hLRtotS->GetNbinsX(),Eff,Pur);
00448   hEffvsPur -> SetName("hEffvsPur");
00449   hEffvsPur -> SetTitle("");
00450   hEffvsPur -> GetXaxis() -> SetTitle((TString)("Efficiency of cut on log combined LR"));
00451   hEffvsPur -> GetYaxis() -> SetTitle((TString)("Purity"));
00452   hEffvsPur -> GetYaxis() -> SetRangeUser(0,1.1);
00453 
00454   hLRValvsPur = new TGraph(hLRtotS->GetNbinsX(),LRVal,Pur);
00455   hLRValvsPur -> SetName("hLRValvsPur");
00456   hLRValvsPur -> SetTitle("");
00457   hLRValvsPur -> GetXaxis() -> SetTitle((TString)("Cut on the log combined LR value"));
00458   hLRValvsPur -> GetYaxis() -> SetTitle((TString)("Purity"));
00459   hLRValvsPur -> GetYaxis() -> SetRangeUser(0,1.1);
00460 
00461   hLRValvsEff = new TGraph(hLRtotS->GetNbinsX(),LRVal,Eff);
00462   hLRValvsEff -> SetName("hLRValvsEff");
00463   hLRValvsEff -> SetTitle("");
00464   hLRValvsEff -> GetXaxis() -> SetTitle((TString)("Cut on the log combined LR value"));
00465   hLRValvsEff -> GetYaxis() -> SetTitle((TString)("Efficiency of cut on log combined LR"));
00466   hLRValvsEff -> GetYaxis() -> SetRangeUser(0,1.1);
00467 
00468 }   
00469 
00470 
00471 
00472 
00473 // member function to calculate the probability for signal given a logLR value
00474 double  LRHelpFunctions::calcProb(double logLR){
00475   return fLRtotSoverSplusB -> Eval(logLR);
00476 }
00477 
00478 
00479 // member to check if a certain S/S+B fit function is read from the root-file
00480 bool    LRHelpFunctions::obsFitIncluded(int o){
00481   bool included = false;
00482   TString obs = "_Obs"; obs += o; obs += "_";
00483   for(size_t f = 0; f<fObsSoverSplusB.size(); f++){    
00484     if(((TString)(fObsSoverSplusB[f]->GetName())).Contains(obs)) included = true;
00485   }
00486   return included;
00487 }
00488 
00489 void LRHelpFunctions::setXlabels(const std::vector<std::string> & xLabels)
00490 {
00491   if (hObsS.size() != xLabels.size()) {
00492     std::cout << "LRHelpFunctions::setXlabels: Number of labels ("
00493          << xLabels.size()<< ") does not match number of obervables("
00494          << hObsS.size()<<").\n";
00495     return;
00496   }
00497   for(size_t i = 0; i<hObsS.size(); ++i){
00498     hObsS[i]->GetXaxis()->SetTitle(TString(xLabels[i]));
00499     hObsB[i]->GetXaxis()->SetTitle(TString(xLabels[i]));
00500     hObsSoverSplusB[i]->GetXaxis()->SetTitle(TString(xLabels[i]));
00501     hObsSoverSplusB[i]->GetYaxis()->SetTitle(TString("S/(S+B)"));
00502   }
00503 }
00504 
00505 void LRHelpFunctions::setYlabels(const std::vector<std::string> & yLabels)
00506 {
00507   if (hObsS.size() != yLabels.size()) {
00508     std::cout << "LRHelpFunctions::setYlabels: Number of labels ("
00509          << yLabels.size()<< ") does not match number of obervables("
00510          << hObsS.size()<<").\n";
00511     return;
00512   }
00513   for(size_t i = 0; i<hObsS.size(); ++i){
00514     hObsS[i]->GetYaxis()->SetTitle(TString(yLabels[i]));
00515     hObsB[i]->GetYaxis()->SetTitle(TString(yLabels[i]));
00516   }
00517 }
00518 
00519 void  LRHelpFunctions::singlePlot(TString fname, int obsNbr, TString extension) {
00520   if (!obsFitIncluded(obsNbr)) return;
00521 
00522   TStyle *tdrStyle = gROOT->GetStyle("tdrStyle");
00523   tdrStyle->SetOptFit(0);
00524   tdrStyle->SetOptStat(0);
00525   tdrStyle->SetOptTitle(0);
00526 //   tdrStyle->SetPadTopMargin(0.01);
00527 //   tdrStyle->SetPadBottomMargin(0.01);
00528 //   tdrStyle->SetPadLeftMargin(0.01);
00529 //   tdrStyle->SetPadRightMargin(0.01);
00530 
00531   TCanvas c2("c2","",600,300);
00532   c2.SetTopMargin(0.01);
00533   c2.SetBottomMargin(0.01);
00534   c2.SetLeftMargin(0.01);
00535   c2.SetRightMargin(0.01);
00536   std::cout <<fname<<std::endl;
00537   c2.Divide(2,1);
00538 
00539   TString obs = "Obs"; obs += obsNbr; obs += "_";
00540   for(size_t o = 0; o<hObsB.size(); ++o){
00541     if(((TString)(hObsB[o]->GetName())).Contains(obs)) {
00542       c2.cd(1);
00543       hObsS[o] -> SetLineColor(2);
00544       if(hObsB[o]->GetMaximum() > hObsS[o]->GetMaximum()){
00545         hObsB[o] -> Draw("hist");
00546         hObsS[o] -> Draw("histsame");
00547       }
00548       else
00549       {
00550         hObsS[o] -> Draw("hist");
00551         hObsB[o] -> Draw("histsame");
00552       }
00553       c2.cd(2);
00554 
00555       hObsSoverSplusB[o] -> Draw();
00556       fObsSoverSplusB[o] -> Draw("same");
00557     }
00558   }
00559   std::cout << fname+"."+extension<<std::endl;
00560   c2.Print(fname+"."+extension);
00561 }
00562 
00563 void  LRHelpFunctions::purityPlot(TString fname, TString extension)
00564 {
00565   TStyle *tdrStyle = gROOT->GetStyle("tdrStyle");
00566   tdrStyle->SetOptFit(0);
00567   tdrStyle->SetOptStat(0);
00568   tdrStyle->SetOptTitle(0);
00569 //   tdrStyle->SetPadTopMargin(0.01);
00570 //   tdrStyle->SetPadBottomMargin(0.01);
00571 //   tdrStyle->SetPadLeftMargin(0.01);
00572 //   tdrStyle->SetPadRightMargin(0.01);
00573 
00574   TCanvas c2("c2","",600,300);
00575   c2.SetTopMargin(0.01);
00576   c2.SetBottomMargin(0.01);
00577   c2.SetLeftMargin(0.01);
00578   c2.SetRightMargin(0.01);
00579   std::cout <<fname<<std::endl;
00580   c2.Divide(2,1);
00581 
00582   c2.cd(1);
00583   hLRValvsPur->Draw("AP");
00584   c2.cd(2);
00585   hEffvsPur->Draw("AP");
00586   c2.Print(fname+"Purity."+extension);
00587 
00588   hLRtotS->GetXaxis()->SetNdivisions(505);
00589   hLRtotB->GetXaxis()->SetNdivisions(505);
00590   TCanvas c3("c2","",300,300);
00591 //   c3.SetTopMargin(0.01);
00592 //   c3.SetBottomMargin(0.01);
00593 //   c3.SetLeftMargin(0.01);
00594 //   c3.SetRightMargin(0.01);
00595   hLRtotS->GetXaxis()->SetTitle("Combined LR");
00596   hLRtotB->GetXaxis()->SetTitle("Combined LR");
00597 
00598     hLRtotB -> Draw();
00599     hLRtotS -> SetLineColor(2);
00600     hLRtotS -> Draw("same");
00601   c3.Print(fname+"Dist."+extension);
00602 
00603 
00604     TCanvas c5("c3","",900,600);
00605     c5.Divide(2,2);
00606     c5.cd(1);
00607     hLRtotB -> Draw();
00608     hLRtotS -> SetLineColor(2);
00609     hLRtotS -> Draw("same");
00610     c5.cd(3);
00611     hLRValvsEff -> Draw("AP");
00612     c5.cd(2);
00613     hEffvsPur -> Draw("AP");
00614     c5.cd(4);
00615     hLRtotSoverSplusB -> Draw();
00616     c5.Print(fname+"all."+extension);
00617 }