CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
LRHelpFunctions Class Reference

Help functionalities to implement and evaluate LR ratio method. More...

#include "TopQuarkAnalysis/TopTools/interface/LRHelpFunctions.h"

Public Member Functions

double calcLRval (const std::vector< double > &)
 
double calcProb (double)
 
double calcPtdrLRval (const std::vector< double > &vals, bool useCorrelation=false)
 
void fillLRBackgroundHist (double, double weight=1.0)
 
void fillLRSignalHist (double, double weight=1.0)
 
void fillToBackgroundHists (const std::vector< double > &, double weight=1.0)
 
void fillToBackgroundHists (int, double, double weight=1.0)
 
void fillToSignalCorrelation (int obsNbr1, double obsVal1, int obsNbr2, double obsVal2, double weight)
 
void fillToSignalHists (const std::vector< double > &, double weight=1.0)
 
void fillToSignalHists (int, double, double weight=1.0)
 
void initLRHistsAndFits (int, double, double, const char *)
 
 LRHelpFunctions ()
 
 LRHelpFunctions (const std::vector< int > &, int, const std::vector< double > &, const std::vector< double > &, const std::vector< const char * > &)
 
 LRHelpFunctions (int, double, double, const char *)
 
void makeAndFitPurityHists ()
 
void makeAndFitSoverSplusBHists ()
 
void normalizeSandBhists ()
 
bool obsFitIncluded (int)
 
void purityPlot (const TString &fname, const TString &extension=TString("eps"))
 
void readObsHistsAndFits (const TString &, const std::vector< int > &, bool)
 
void recreateFitFct (const std::vector< int > &obsNr, const std::vector< const char * > &functions)
 
void setObsFitParameters (int obs, const std::vector< double > &)
 
void setXlabels (const std::vector< std::string > &xLabels)
 
void setYlabels (const std::vector< std::string > &yLabels)
 
void singlePlot (const TString &fname, int obsNbr, const TString &extension=TString("eps"))
 
void storeControlPlots (const TString &)
 
void storeToROOTfile (const TString &)
 
 ~LRHelpFunctions ()
 

Private Attributes

bool constructPurity
 
TF1 * fLRtotSoverSplusB
 
std::vector< TF1 * > fObsSoverSplusB
 
TGraph * hEffvsPur
 
TH1F * hLRtotB
 
TH1F * hLRtotS
 
TH1F * hLRtotSoverSplusB
 
TGraph * hLRValvsEff
 
TGraph * hLRValvsPur
 
std::vector< TH1F * > hObsB
 
std::vector< TH2F * > hObsCorr
 
std::vector< TH1F * > hObsS
 
std::vector< TH1F * > hObsSoverSplusB
 
std::vector< int > obsNumbers
 

Detailed Description

Help functionalities to implement and evaluate LR ratio method.

Author
Jan Heyninck
Version
Id:
LRHelpFunctions.h,v 1.11 2013/05/30 20:47:31 gartung Exp

Definition at line 37 of file LRHelpFunctions.h.

Constructor & Destructor Documentation

LRHelpFunctions::LRHelpFunctions ( )

Definition at line 5 of file LRHelpFunctions.cc.

References constructPurity, and plotscripts::setTDRStyle().

5  {
6  constructPurity = false;
7  setTDRStyle();
8  gStyle->SetCanvasDefW(900);
9 }
def setTDRStyle
Definition: plotscripts.py:87
LRHelpFunctions::LRHelpFunctions ( const std::vector< int > &  obsNr,
int  nrBins,
const std::vector< double > &  obsMin,
const std::vector< double > &  obsMax,
const std::vector< const char * > &  functions 
)

Definition at line 12 of file LRHelpFunctions.cc.

References constructPurity, fObsSoverSplusB, hObsB, hObsCorr, hObsS, hObsSoverSplusB, python.connectstrParser::o, obsNumbers, and plotscripts::setTDRStyle().

12  {
13  obsNumbers = obsNr;
14  constructPurity = false;
15  setTDRStyle();
16  gStyle->SetCanvasDefW(900);
17  for(size_t o=0; o<obsNr.size(); o++){
18  // create Signal, background and s/(s+b) histogram
19  TString htS = "Obs"; htS += obsNr[o]; htS += "_S";
20  TString htB = "Obs"; htB += obsNr[o]; htB += "_B";
21  TString htSB = "Obs"; htSB += obsNr[o]; htSB += "_SoverSplusB";
22  hObsS.push_back( new TH1F(htS,htS,nrBins,obsMin[o],obsMax[o]) );
23  hObsB.push_back( new TH1F(htB,htB,nrBins,obsMin[o],obsMax[o]) );
24  hObsSoverSplusB.push_back( new TH1F(htSB,htSB,nrBins,obsMin[o],obsMax[o]) );
25 
26  // create the correlation 2D plots among the observables (for signal)
27  for (size_t o2 = o+1; o2 < obsNr.size(); o2++) {
28  TString hcorr = "Corr_Obs"; hcorr += obsNr[o]; hcorr += "_Obs"; hcorr += obsNr[o2];
29  hObsCorr.push_back( new TH2F(hcorr,hcorr,nrBins,obsMin[o],obsMax[o],nrBins,obsMin[o2],obsMax[o2]) );
30  }
31 
32  // create fit functions
33  TString ftSB = "F_Obs"; ftSB += obsNr[o]; ftSB += "_SoverSplusB";
34  fObsSoverSplusB.push_back( new TF1(ftSB,functions[o],hObsS[o]->GetXaxis()->GetXmin(),hObsS[o]->GetXaxis()->GetXmax()) );
35  }
36 }
std::vector< TH1F * > hObsSoverSplusB
def setTDRStyle
Definition: plotscripts.py:87
std::vector< TH1F * > hObsS
std::vector< TH1F * > hObsB
std::vector< TH2F * > hObsCorr
std::vector< int > obsNumbers
std::vector< TF1 * > fObsSoverSplusB
LRHelpFunctions::LRHelpFunctions ( int  nrLRbins,
double  LRmin,
double  LRmax,
const char *  LRfunction 
)

Definition at line 47 of file LRHelpFunctions.cc.

References constructPurity, fLRtotSoverSplusB, hLRtotB, hLRtotS, hLRtotSoverSplusB, and plotscripts::setTDRStyle().

47  {
48  constructPurity = true;
49  setTDRStyle();
50  gStyle->SetCanvasDefW(900);
51 
52  // create LR histograms
53  hLRtotS = new TH1F("hLRtotS","hLRtotS",nrLRbins,LRmin,LRmax);
54  hLRtotS->GetXaxis()->SetTitle("Combined LR");
55  hLRtotB = new TH1F("hLRtotB","hLRtotB",nrLRbins,LRmin,LRmax);
56  hLRtotB->GetXaxis()->SetTitle("Combined LR");
57  hLRtotSoverSplusB = new TH1F("hLRtotSoverSplusB","hLRtotSoverSplusB",nrLRbins,LRmin,LRmax);
58 
59  // create LR fit function
60  fLRtotSoverSplusB = new TF1("fLRtotSoverSplusB",LRfunction,LRmin,LRmax);
61 }
def setTDRStyle
Definition: plotscripts.py:87
TH1F * hLRtotSoverSplusB
LRHelpFunctions::~LRHelpFunctions ( )

Definition at line 65 of file LRHelpFunctions.cc.

65 {}

Member Function Documentation

double LRHelpFunctions::calcLRval ( const std::vector< double > &  vals)

Definition at line 361 of file LRHelpFunctions.cc.

References fObsSoverSplusB, create_public_lumi_plots::log, and python.connectstrParser::o.

Referenced by TtHadLRJetCombCalc::operator()(), TtHadLRSignalSelCalc::operator()(), TtSemiLRJetCombCalc::operator()(), and TtSemiLRSignalSelCalc::operator()().

361  {
362  double logLR = 0.;
363  for(size_t o=0; o<fObsSoverSplusB.size(); o++){
364  double SoverSplusN = fObsSoverSplusB[o]->Eval(vals[o]);
365  double SoverN = 1./((1./SoverSplusN) - 1.);
366  logLR += log(SoverN);
367  }
368  return logLR;
369 }
std::vector< TF1 * > fObsSoverSplusB
double LRHelpFunctions::calcProb ( double  logLR)
double LRHelpFunctions::calcPtdrLRval ( const std::vector< double > &  vals,
bool  useCorrelation = false 
)

Definition at line 375 of file LRHelpFunctions.cc.

References corr, fObsSoverSplusB, hObsCorr, i, j, max(), min, python.connectstrParser::o, obsNumbers, and funct::pow().

375  {
376  double logLR = 1.;
377  for(size_t o=0; o<fObsSoverSplusB.size(); o++){
378  double SoverSplusN = fObsSoverSplusB[o]->Eval(vals[o]);
379  if (SoverSplusN<0.0001) SoverSplusN=0.0001;
380  if (useCorrelation){
381  double corr = 0;
382  for(size_t j=0; j<fObsSoverSplusB.size(); j++){
383  if (o==j) ++corr;
384  else {
385  TString hcorr = "Corr_Obs"; hcorr += std::min(obsNumbers[o],obsNumbers[j]) ; hcorr += "_Obs"; hcorr += std::max(obsNumbers[o],obsNumbers[j]);
386  for(size_t i=0; i<hObsCorr.size(); i++) {
387  if(((TString)(hObsCorr[i]->GetName())) == hcorr)
388  corr += fabs(hObsCorr[i]->GetCorrelationFactor());
389  }
390  }
391  }
392  logLR *= pow(SoverSplusN/fObsSoverSplusB[o]->GetMaximum(), 1./corr);
393  }
394  else logLR *= SoverSplusN/fObsSoverSplusB[o]->GetMaximum();
395  }
396  //std::cout <<logLR<<std::endl;
397  return logLR;
398 }
int i
Definition: DBlmapReader.cc:9
#define min(a, b)
Definition: mlp_lapack.h:161
std::vector< TH2F * > hObsCorr
const T & max(const T &a, const T &b)
int j
Definition: DBlmapReader.cc:9
JetCorrectorParameters corr
Definition: classes.h:11
std::vector< int > obsNumbers
std::vector< TF1 * > fObsSoverSplusB
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
void LRHelpFunctions::fillLRBackgroundHist ( double  val,
double  weight = 1.0 
)

Definition at line 410 of file LRHelpFunctions.cc.

References HcalObjRepresent::Fill(), and hLRtotB.

410  {
411  // std::cout<< "@@@===> LRHelpf Backgroud LRval = "<< val<<std::endl;
412  hLRtotB -> Fill(val, weight);
413 }
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
int weight
Definition: histoStyle.py:50
void LRHelpFunctions::fillLRSignalHist ( double  val,
double  weight = 1.0 
)

Definition at line 402 of file LRHelpFunctions.cc.

References HcalObjRepresent::Fill(), and hLRtotS.

402  {
403  // std::cout<< "@@@===> LRHelpf Signal LRval = "<< val<<std::endl;
404  hLRtotS -> Fill(val, weight);
405 }
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
int weight
Definition: histoStyle.py:50
void LRHelpFunctions::fillToBackgroundHists ( const std::vector< double > &  obsVals,
double  weight = 1.0 
)

Definition at line 138 of file LRHelpFunctions.cc.

References HcalObjRepresent::Fill(), hObsB, and python.connectstrParser::o.

138  {
139  for(size_t o=0; o<obsVals.size(); o++) hObsB[o]->Fill(obsVals[o], weight);
140 }
std::vector< TH1F * > hObsB
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
int weight
Definition: histoStyle.py:50
void LRHelpFunctions::fillToBackgroundHists ( int  obsNbr,
double  obsVal,
double  weight = 1.0 
)

Definition at line 143 of file LRHelpFunctions.cc.

References f, and hObsB.

143  {
144  TString obs = "Obs"; obs += obsNbr; obs += "_";
145  for(size_t f = 0; f<hObsB.size(); f++){
146  if(((TString)(hObsB[f]->GetName())).Contains(obs)) {
147  hObsB[f]->Fill(obsVal, weight);
148  return;
149  }
150  }
151 }
std::vector< TH1F * > hObsB
double f[11][100]
int weight
Definition: histoStyle.py:50
void LRHelpFunctions::fillToSignalCorrelation ( int  obsNbr1,
double  obsVal1,
int  obsNbr2,
double  obsVal2,
double  weight 
)

Definition at line 121 of file LRHelpFunctions.cc.

References HcalObjRepresent::Fill(), hObsCorr, i, max(), and min.

122  {
123  TString hcorr = "Corr_Obs"; hcorr += std::min(obsNbr1,obsNbr2) ; hcorr += "_Obs"; hcorr += std::max(obsNbr1, obsNbr2);
124  for(size_t i=0; i<hObsCorr.size(); i++) {
125  if(((TString)(hObsCorr[i]->GetName())) == hcorr) {
126  if (obsNbr1 < obsNbr2) {
127  hObsCorr[i] -> Fill(obsVal1,obsVal2, weight);
128  } else {
129  hObsCorr[i] -> Fill(obsVal2,obsVal1, weight);
130  }
131  return;
132  }
133  }
134 }
int i
Definition: DBlmapReader.cc:9
#define min(a, b)
Definition: mlp_lapack.h:161
std::vector< TH2F * > hObsCorr
const T & max(const T &a, const T &b)
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
int weight
Definition: histoStyle.py:50
void LRHelpFunctions::fillToSignalHists ( const std::vector< double > &  obsVals,
double  weight = 1.0 
)

Definition at line 98 of file LRHelpFunctions.cc.

References HcalObjRepresent::Fill(), hObsCorr, hObsS, and python.connectstrParser::o.

98  {
99  int hIndex = 0;
100  for(size_t o=0; o<obsVals.size(); o++) {
101  hObsS[o]->Fill(obsVals[o], weight);
102  for(size_t o2=o+1; o2<obsVals.size(); o2++) {
103  hObsCorr[hIndex] -> Fill(obsVals[o],obsVals[o2], weight);
104  ++hIndex;
105  }
106  }
107 }
std::vector< TH1F * > hObsS
std::vector< TH2F * > hObsCorr
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
int weight
Definition: histoStyle.py:50
void LRHelpFunctions::fillToSignalHists ( int  obsNbr,
double  obsVal,
double  weight = 1.0 
)

Definition at line 110 of file LRHelpFunctions.cc.

References f, and hObsS.

110  {
111  TString obs = "Obs"; obs += obsNbr; obs += "_";
112  for(size_t f = 0; f<hObsS.size(); f++){
113  if(((TString)(hObsS[f]->GetName())).Contains(obs)) {
114  hObsS[f]->Fill(obsVal, weight);
115  return;
116  }
117  }
118 }
std::vector< TH1F * > hObsS
double f[11][100]
int weight
Definition: histoStyle.py:50
void LRHelpFunctions::initLRHistsAndFits ( int  nrLRbins,
double  LRmin,
double  LRmax,
const char *  LRfunction 
)

Definition at line 69 of file LRHelpFunctions.cc.

References constructPurity, fLRtotSoverSplusB, hLRtotB, hLRtotS, and hLRtotSoverSplusB.

69  {
70  constructPurity = true;
71  // create LR histograms
72  hLRtotS = new TH1F("hLRtotS","hLRtotS",nrLRbins,LRmin,LRmax);
73  hLRtotB = new TH1F("hLRtotB","hLRtotB",nrLRbins,LRmin,LRmax);
74  hLRtotSoverSplusB = new TH1F("hLRtotSoverSplusB","hLRtotSoverSplusB",nrLRbins,LRmin,LRmax);
75  // create LR fit function
76  fLRtotSoverSplusB = new TF1("fLRtotSoverSplusB",LRfunction,LRmin,LRmax);
77 }
TH1F * hLRtotSoverSplusB
void LRHelpFunctions::makeAndFitPurityHists ( )

Definition at line 418 of file LRHelpFunctions.cc.

References b, gather_cfg::cout, GOODCOLL_filter_cfg::cut, fLRtotSoverSplusB, hEffvsPur, hLRtotB, hLRtotS, hLRtotSoverSplusB, hLRValvsEff, hLRValvsPur, funct::pow(), and mathSSE::sqrt().

418  {
419 
420  for(int b=0; b<=hLRtotS->GetNbinsX(); b++) {
421  float Sint = hLRtotS->Integral(b, hLRtotS->GetNbinsX()+1);
422  float Bint = hLRtotB->Integral(b, hLRtotB->GetNbinsX()+1);
423  if (Sint + Bint > 0) {
424  hLRtotSoverSplusB->SetBinContent(b,1. * Sint / (Sint + Bint));
425  hLRtotSoverSplusB->SetBinError(b,sqrt((pow(Sint*sqrt(Bint),2)+pow(Bint*sqrt(Sint),2)))/pow((Sint+Bint),2));
426  }
427  }
428 
429  hLRtotS->GetXaxis()->SetTitle("Combined LR ratio");
430  hLRtotB->GetXaxis()->SetTitle("Combined LR ratio");
431  hLRtotSoverSplusB->GetXaxis()->SetTitle("Cut on Combined LR");
432  hLRtotSoverSplusB->GetYaxis()->SetTitle("Purity");
433 
434  hLRtotSoverSplusB -> Fit(fLRtotSoverSplusB->GetName(),"RQ");
435  double totSignal = hLRtotS->Integral(0,hLRtotS->GetNbinsX()+1);
436  double Eff[200], Pur[200], LRVal[200];
437  if (hLRtotS->GetNbinsX()>200) {
438  std::cout << "Number of bins of LR histograms can not execeed 200!"<<std::endl;
439  return;
440  }
441  for(int cut=0; (cut<=hLRtotS->GetNbinsX())&&(cut<200) ; cut ++){
442  double LRcutVal = hLRtotS->GetBinLowEdge(cut);
443  Eff[cut] = hLRtotS->Integral(cut,hLRtotS->GetNbinsX()+1)/totSignal;
444  Pur[cut] = fLRtotSoverSplusB->Eval(LRcutVal);
445  LRVal[cut] = LRcutVal;
446  }
447  hEffvsPur = new TGraph(hLRtotS->GetNbinsX(),Eff,Pur);
448  hEffvsPur -> SetName("hEffvsPur");
449  hEffvsPur -> SetTitle("");
450  hEffvsPur -> GetXaxis() -> SetTitle((TString)("Efficiency of cut on log combined LR"));
451  hEffvsPur -> GetYaxis() -> SetTitle((TString)("Purity"));
452  hEffvsPur -> GetYaxis() -> SetRangeUser(0,1.1);
453 
454  hLRValvsPur = new TGraph(hLRtotS->GetNbinsX(),LRVal,Pur);
455  hLRValvsPur -> SetName("hLRValvsPur");
456  hLRValvsPur -> SetTitle("");
457  hLRValvsPur -> GetXaxis() -> SetTitle((TString)("Cut on the log combined LR value"));
458  hLRValvsPur -> GetYaxis() -> SetTitle((TString)("Purity"));
459  hLRValvsPur -> GetYaxis() -> SetRangeUser(0,1.1);
460 
461  hLRValvsEff = new TGraph(hLRtotS->GetNbinsX(),LRVal,Eff);
462  hLRValvsEff -> SetName("hLRValvsEff");
463  hLRValvsEff -> SetTitle("");
464  hLRValvsEff -> GetXaxis() -> SetTitle((TString)("Cut on the log combined LR value"));
465  hLRValvsEff -> GetYaxis() -> SetTitle((TString)("Efficiency of cut on log combined LR"));
466  hLRValvsEff -> GetYaxis() -> SetRangeUser(0,1.1);
467 
468 }
TGraph * hLRValvsEff
TGraph * hLRValvsPur
Definition: Fit.h:34
T sqrt(T t)
Definition: SSEVec.h:48
double b
Definition: hdecay.h:120
tuple cout
Definition: gather_cfg.py:121
TH1F * hLRtotSoverSplusB
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
void LRHelpFunctions::makeAndFitSoverSplusBHists ( )

Definition at line 176 of file LRHelpFunctions.cc.

References b, error, fObsSoverSplusB, hObsB, hObsS, hObsSoverSplusB, python.connectstrParser::o, funct::pow(), and mathSSE::sqrt().

176  {
177  for(size_t o=0; o<hObsS.size(); o++){
178  for(int b=0; b <= hObsS[o]->GetNbinsX()+1; b++){
179  if ((hObsS[o]->GetBinContent(b)+ hObsB[o]->GetBinContent(b)) > 0) {
180  hObsSoverSplusB[o]->SetBinContent(b, hObsS[o]->GetBinContent(b) / (hObsS[o]->GetBinContent(b) + hObsB[o]->GetBinContent(b)));
181  double error = sqrt( pow(hObsS[o]->GetBinError(b) * hObsB[o]->GetBinContent(b)/pow(hObsS[o]->GetBinContent(b) + hObsB[o]->GetBinContent(b),2),2)
182  + pow(hObsB[o]->GetBinError(b) * hObsS[o]->GetBinContent(b)/pow(hObsS[o]->GetBinContent(b) + hObsB[o]->GetBinContent(b),2),2) );
183  hObsSoverSplusB[o]->SetBinError(b,error);
184  }
185  }
186  hObsSoverSplusB[o]->Fit(fObsSoverSplusB[o]->GetName(),"RQ");
187  }
188 }
std::vector< TH1F * > hObsSoverSplusB
std::vector< TH1F * > hObsS
std::vector< TH1F * > hObsB
T sqrt(T t)
Definition: SSEVec.h:48
double b
Definition: hdecay.h:120
std::vector< TF1 * > fObsSoverSplusB
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
void LRHelpFunctions::normalizeSandBhists ( )

Definition at line 154 of file LRHelpFunctions.cc.

References b, hObsB, hObsS, i, and python.connectstrParser::o.

154  {
155  for(size_t o=0; o<hObsS.size(); o++){
156  // count entries in each histo. Do it this way instead of GetEntries method,
157  // since the latter does not account for the weights.
158  double nrSignEntries = 0, nrBackEntries = 0;
159  for (int i=0; i<=hObsS[o]->GetNbinsX()+1; i++){
160  nrSignEntries += hObsS[o]->GetBinContent(i);
161  nrBackEntries += hObsB[o]->GetBinContent(i);
162  }
163  for(int b=0; b <= hObsS[o]->GetNbinsX()+1; b++){
164  hObsS[o]->SetBinContent(b,hObsS[o]->GetBinContent(b)/(nrSignEntries));
165  hObsB[o]->SetBinContent(b,hObsB[o]->GetBinContent(b)/(nrBackEntries));
166  hObsS[o]->SetBinError(b,hObsS[o]->GetBinError(b)/(nrSignEntries));
167  hObsB[o]->SetBinError(b,hObsB[o]->GetBinError(b)/(nrBackEntries));
168  }
169  //std::cout<<"Integral for obs"<<o<<" S: "<<hObsS[o]->Integral(0,10000)<<" & obs"<<o<<" B: "<<hObsB[o]->Integral(0,10000)<<std::endl;
170  }
171 }
int i
Definition: DBlmapReader.cc:9
std::vector< TH1F * > hObsS
std::vector< TH1F * > hObsB
double b
Definition: hdecay.h:120
bool LRHelpFunctions::obsFitIncluded ( int  o)

Definition at line 480 of file LRHelpFunctions.cc.

References f, fObsSoverSplusB, and python.connectstrParser::o.

Referenced by TtHadLRJetCombCalc::operator()(), TtHadLRSignalSelCalc::operator()(), TtSemiLRJetCombCalc::operator()(), TtSemiLRSignalSelCalc::operator()(), and singlePlot().

480  {
481  bool included = false;
482  TString obs = "_Obs"; obs += o; obs += "_";
483  for(size_t f = 0; f<fObsSoverSplusB.size(); f++){
484  if(((TString)(fObsSoverSplusB[f]->GetName())).Contains(obs)) included = true;
485  }
486  return included;
487 }
double f[11][100]
std::vector< TF1 * > fObsSoverSplusB
void LRHelpFunctions::purityPlot ( const TString &  fname,
const TString &  extension = TString("eps") 
)

Definition at line 563 of file LRHelpFunctions.cc.

References gather_cfg::cout, hEffvsPur, hLRtotB, hLRtotS, hLRtotSoverSplusB, hLRValvsEff, hLRValvsPur, and plotscripts::tdrStyle.

564 {
565  TStyle *tdrStyle = gROOT->GetStyle("tdrStyle");
566  tdrStyle->SetOptFit(0);
567  tdrStyle->SetOptStat(0);
568  tdrStyle->SetOptTitle(0);
569 // tdrStyle->SetPadTopMargin(0.01);
570 // tdrStyle->SetPadBottomMargin(0.01);
571 // tdrStyle->SetPadLeftMargin(0.01);
572 // tdrStyle->SetPadRightMargin(0.01);
573 
574  TCanvas c2("c2","",600,300);
575  c2.SetTopMargin(0.01);
576  c2.SetBottomMargin(0.01);
577  c2.SetLeftMargin(0.01);
578  c2.SetRightMargin(0.01);
579  std::cout <<fname<<std::endl;
580  c2.Divide(2,1);
581 
582  c2.cd(1);
583  hLRValvsPur->Draw("AP");
584  c2.cd(2);
585  hEffvsPur->Draw("AP");
586  c2.Print(fname+"Purity."+extension);
587 
588  hLRtotS->GetXaxis()->SetNdivisions(505);
589  hLRtotB->GetXaxis()->SetNdivisions(505);
590  TCanvas c3("c2","",300,300);
591 // c3.SetTopMargin(0.01);
592 // c3.SetBottomMargin(0.01);
593 // c3.SetLeftMargin(0.01);
594 // c3.SetRightMargin(0.01);
595  hLRtotS->GetXaxis()->SetTitle("Combined LR");
596  hLRtotB->GetXaxis()->SetTitle("Combined LR");
597 
598  hLRtotB -> Draw();
599  hLRtotS -> SetLineColor(2);
600  hLRtotS -> Draw("same");
601  c3.Print(fname+"Dist."+extension);
602 
603 
604  TCanvas c5("c3","",900,600);
605  c5.Divide(2,2);
606  c5.cd(1);
607  hLRtotB -> Draw();
608  hLRtotS -> SetLineColor(2);
609  hLRtotS -> Draw("same");
610  c5.cd(3);
611  hLRValvsEff -> Draw("AP");
612  c5.cd(2);
613  hEffvsPur -> Draw("AP");
614  c5.cd(4);
615  hLRtotSoverSplusB -> Draw();
616  c5.Print(fname+"all."+extension);
617 }
TGraph * hLRValvsEff
TGraph * hLRValvsPur
string fname
main script
tuple cout
Definition: gather_cfg.py:121
TH1F * hLRtotSoverSplusB
void LRHelpFunctions::readObsHistsAndFits ( const TString &  fileName,
const std::vector< int > &  observables,
bool  readLRplots 
)

Definition at line 192 of file LRHelpFunctions.cc.

References constructPurity, gather_cfg::cout, fLRtotSoverSplusB, fObsSoverSplusB, hLRtotB, hLRtotS, hLRtotSoverSplusB, hObsB, hObsCorr, hObsS, hObsSoverSplusB, list(), and obsNumbers.

192  {
193  hObsS.clear();
194  hObsB.clear();
195  hObsSoverSplusB.clear();
196  fObsSoverSplusB.clear();
197  TFile *fitFile = new TFile(fileName, "READ");
198  if(observables[0] == -1){
199  std::cout<<" ... will read hists and fit for all available observables in file "<<fileName<<std::endl;
200  TList *list = fitFile->GetListOfKeys();
201  TIter next(list);
202  TKey *el;
203  while ((el = (TKey*)next())) {
204  TString keyName = el->GetName();
205  if(keyName.Contains("F_") && keyName.Contains("_SoverSplusB")){
206  fObsSoverSplusB.push_back( new TF1(*((TF1*) el -> ReadObj())));
207  }
208  else if(keyName.Contains("_SoverSplusB")){
209  hObsSoverSplusB.push_back( new TH1F(*((TH1F*) el -> ReadObj())));
210  }
211  else if(keyName.Contains("_S")){
212  hObsS.push_back( new TH1F(*((TH1F*) el -> ReadObj())));
213  }
214  else if(keyName.Contains("_B")){
215  hObsB.push_back( new TH1F(*((TH1F*) el -> ReadObj())));
216  }
217  else if(keyName.Contains("Corr")){
218  hObsCorr.push_back( new TH2F(*((TH2F*) el -> ReadObj())));
219  }
220  }
221  }
222  else{
223  obsNumbers = observables;
224  for(unsigned int obs = 0; obs < observables.size(); obs++){
225  std::cout<<" ... will read hists and fit for obs "<<observables[obs]<<" from file "<<fileName<<std::endl;
226  TString hStitle = "Obs"; hStitle += observables[obs]; hStitle += "_S";
227  hObsS.push_back( new TH1F(*((TH1F*)fitFile->GetKey(hStitle)->ReadObj())));
228  TString hBtitle = "Obs"; hBtitle += observables[obs]; hBtitle += "_B";
229  hObsB.push_back( new TH1F(*((TH1F*)fitFile->GetKey(hBtitle)->ReadObj())));
230  TString hSBtitle = "Obs"; hSBtitle += observables[obs]; hSBtitle += "_SoverSplusB";
231  TString fSBtitle = "F_"; fSBtitle += hSBtitle;
232  hObsSoverSplusB.push_back( new TH1F(*((TH1F*)fitFile->GetKey(hSBtitle)->ReadObj())));
233  fObsSoverSplusB.push_back( new TF1(*((TF1*)fitFile->GetKey(fSBtitle)->ReadObj())));
234  for(unsigned int obs2 = obs+1; obs2 < observables.size(); obs2++){
235  TString hCorrtitle = "Corr_Obs"; hCorrtitle += observables[obs];
236  hCorrtitle += "_Obs"; hCorrtitle += observables[obs2];
237  hObsCorr.push_back( new TH2F(*((TH2F*)fitFile->GetKey(hCorrtitle)->ReadObj())));
238  }
239  }
240  }
241 
242  if(readLRplots){
243  constructPurity = true;
244  std::cout<<" ... will LR s and B histos from file "<<fileName<<std::endl;
245  hLRtotS = new TH1F(*((TH1F*)fitFile->GetKey("hLRtotS")->ReadObj()));
246  hLRtotB = new TH1F(*((TH1F*)fitFile->GetKey("hLRtotB")->ReadObj()));
247  hLRtotSoverSplusB = new TH1F(*((TH1F*)fitFile->GetKey("hLRtotSoverSplusB")->ReadObj()));
248  fLRtotSoverSplusB = new TF1(*((TF1*)fitFile->GetKey("fLRtotSoverSplusB") -> ReadObj()));
249  }
250 }
std::vector< TH1F * > hObsSoverSplusB
std::vector< TH1F * > hObsS
std::vector< TH1F * > hObsB
std::vector< TH2F * > hObsCorr
std::vector< int > obsNumbers
std::vector< TF1 * > fObsSoverSplusB
tuple cout
Definition: gather_cfg.py:121
TH1F * hLRtotSoverSplusB
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger list("!*","!HLTx*"if it matches 2 triggers or more) will accept the event if all the matching triggers are FAIL.It will reject the event if any of the triggers are PASS or EXCEPTION(this matches the behavior of"!*"before the partial wildcard feature was incorporated).Triggers which are in the READY state are completely ignored.(READY should never be returned since the trigger paths have been run
void LRHelpFunctions::recreateFitFct ( const std::vector< int > &  obsNr,
const std::vector< const char * > &  functions 
)

Definition at line 38 of file LRHelpFunctions.cc.

References fObsSoverSplusB, hObsS, and python.connectstrParser::o.

38  {
39  if (!fObsSoverSplusB.empty()) fObsSoverSplusB.clear();
40  for(size_t o=0; o<obsNr.size(); o++){
41  // create fit functions
42  TString ftSB = "F_Obs"; ftSB += obsNr[o]; ftSB += "_SoverSplusB";
43  fObsSoverSplusB.push_back( new TF1(ftSB,functions[o],hObsS[o]->GetXaxis()->GetXmin(),hObsS[o]->GetXaxis()->GetXmax()) );
44  }
45 }
std::vector< TH1F * > hObsS
std::vector< TF1 * > fObsSoverSplusB
void LRHelpFunctions::setObsFitParameters ( int  obs,
const std::vector< double > &  fitPars 
)

Definition at line 82 of file LRHelpFunctions.cc.

References fObsSoverSplusB, and AlCaHLTBitMon_ParallelJobs::p.

82  {
83  for(size_t fit=0; fit<fObsSoverSplusB.size(); fit++){
84  TString fn = "_Obs"; fn += obs;
85  if(((TString)fObsSoverSplusB[fit]->GetName()).Contains(fn)){
86  for(size_t p=0; p<fitPars.size(); p++){
87  fObsSoverSplusB[fit]->SetParameter(p,fitPars[p]);
88  }
89  }
90  }
91 }
std::vector< TF1 * > fObsSoverSplusB
void LRHelpFunctions::setXlabels ( const std::vector< std::string > &  xLabels)

Definition at line 489 of file LRHelpFunctions.cc.

References gather_cfg::cout, hObsB, hObsS, hObsSoverSplusB, and i.

490 {
491  if (hObsS.size() != xLabels.size()) {
492  std::cout << "LRHelpFunctions::setXlabels: Number of labels ("
493  << xLabels.size()<< ") does not match number of obervables("
494  << hObsS.size()<<").\n";
495  return;
496  }
497  for(size_t i = 0; i<hObsS.size(); ++i){
498  hObsS[i]->GetXaxis()->SetTitle(TString(xLabels[i]));
499  hObsB[i]->GetXaxis()->SetTitle(TString(xLabels[i]));
500  hObsSoverSplusB[i]->GetXaxis()->SetTitle(TString(xLabels[i]));
501  hObsSoverSplusB[i]->GetYaxis()->SetTitle(TString("S/(S+B)"));
502  }
503 }
std::vector< TH1F * > hObsSoverSplusB
int i
Definition: DBlmapReader.cc:9
std::vector< TH1F * > hObsS
std::vector< TH1F * > hObsB
tuple cout
Definition: gather_cfg.py:121
void LRHelpFunctions::setYlabels ( const std::vector< std::string > &  yLabels)

Definition at line 505 of file LRHelpFunctions.cc.

References gather_cfg::cout, hObsB, hObsS, and i.

506 {
507  if (hObsS.size() != yLabels.size()) {
508  std::cout << "LRHelpFunctions::setYlabels: Number of labels ("
509  << yLabels.size()<< ") does not match number of obervables("
510  << hObsS.size()<<").\n";
511  return;
512  }
513  for(size_t i = 0; i<hObsS.size(); ++i){
514  hObsS[i]->GetYaxis()->SetTitle(TString(yLabels[i]));
515  hObsB[i]->GetYaxis()->SetTitle(TString(yLabels[i]));
516  }
517 }
int i
Definition: DBlmapReader.cc:9
std::vector< TH1F * > hObsS
std::vector< TH1F * > hObsB
tuple cout
Definition: gather_cfg.py:121
void LRHelpFunctions::singlePlot ( const TString &  fname,
int  obsNbr,
const TString &  extension = TString("eps") 
)

Definition at line 519 of file LRHelpFunctions.cc.

References gather_cfg::cout, fObsSoverSplusB, hObsB, hObsS, hObsSoverSplusB, python.connectstrParser::o, obsFitIncluded(), and plotscripts::tdrStyle.

519  {
520  if (!obsFitIncluded(obsNbr)) return;
521 
522  TStyle *tdrStyle = gROOT->GetStyle("tdrStyle");
523  tdrStyle->SetOptFit(0);
524  tdrStyle->SetOptStat(0);
525  tdrStyle->SetOptTitle(0);
526 // tdrStyle->SetPadTopMargin(0.01);
527 // tdrStyle->SetPadBottomMargin(0.01);
528 // tdrStyle->SetPadLeftMargin(0.01);
529 // tdrStyle->SetPadRightMargin(0.01);
530 
531  TCanvas c2("c2","",600,300);
532  c2.SetTopMargin(0.01);
533  c2.SetBottomMargin(0.01);
534  c2.SetLeftMargin(0.01);
535  c2.SetRightMargin(0.01);
536  std::cout <<fname<<std::endl;
537  c2.Divide(2,1);
538 
539  TString obs = "Obs"; obs += obsNbr; obs += "_";
540  for(size_t o = 0; o<hObsB.size(); ++o){
541  if(((TString)(hObsB[o]->GetName())).Contains(obs)) {
542  c2.cd(1);
543  hObsS[o] -> SetLineColor(2);
544  if(hObsB[o]->GetMaximum() > hObsS[o]->GetMaximum()){
545  hObsB[o] -> Draw("hist");
546  hObsS[o] -> Draw("histsame");
547  }
548  else
549  {
550  hObsS[o] -> Draw("hist");
551  hObsB[o] -> Draw("histsame");
552  }
553  c2.cd(2);
554 
555  hObsSoverSplusB[o] -> Draw();
556  fObsSoverSplusB[o] -> Draw("same");
557  }
558  }
559  std::cout << fname+"."+extension<<std::endl;
560  c2.Print(fname+"."+extension);
561 }
std::vector< TH1F * > hObsSoverSplusB
std::vector< TH1F * > hObsS
std::vector< TH1F * > hObsB
bool obsFitIncluded(int)
string fname
main script
std::vector< TF1 * > fObsSoverSplusB
tuple cout
Definition: gather_cfg.py:121
void LRHelpFunctions::storeControlPlots ( const TString &  fname)

Definition at line 291 of file LRHelpFunctions.cc.

References trackerHits::c, constructPurity, fObsSoverSplusB, hEffvsPur, hLRtotB, hLRtotS, hLRtotSoverSplusB, hLRValvsEff, hLRValvsPur, hObsB, hObsCorr, hObsS, hObsSoverSplusB, and python.connectstrParser::o.

291  {
292  TCanvas c("dummy","",1);
293  c.Print(fname + "[","landscape");
294  for(size_t o=0; o<hObsS.size(); o++) {
295  TCanvas c2("c2","",1);
296  c2.Divide(2,1);
297  c2.cd(1);
298  hObsS[o] -> SetLineColor(2);
299  if(hObsB[o]->GetMaximum() > hObsS[o]->GetMaximum()){
300  hObsB[o] -> Draw("hist");
301  hObsS[o] -> Draw("histsame");
302  }
303  else
304  {
305  hObsS[o] -> Draw("hist");
306  hObsB[o] -> Draw("histsame");
307  }
308  c2.cd(2);
309  hObsSoverSplusB[o] -> Draw();
310  fObsSoverSplusB[o] -> Draw("same");
311  c2.Print(fname,"Landscape");
312  }
313 
314  int hIndex = 0;
315  for(size_t o=0; o<hObsS.size(); o++) {
316  for(size_t o2=o+1; o2<hObsS.size(); o2++) {
317  TCanvas cc("cc","",1);
318  hObsCorr[hIndex] -> Draw("box");
319  TPaveText pt(0.5,0.87,0.98,0.93,"blNDC");
320  pt.SetFillStyle(1);
321  pt.SetFillColor(0);
322  pt.SetBorderSize(0);
323  TString tcorr = "Corr. of "; tcorr += (int)(100.*hObsCorr[hIndex] -> GetCorrelationFactor()); tcorr += " %";
324  //TText *text = pt.AddText(tcorr);
325  pt.Draw("same");
326  ++hIndex;
327  cc.Print(fname,"Landscape");
328  }
329  }
330 
331  if(constructPurity){
332  TCanvas c3("c3","",1);
333  c3.Divide(2,1);
334  c3.cd(1);
335  hLRtotB -> Draw();
336  hLRtotS -> SetLineColor(2);
337  hLRtotS -> Draw("same");
338  c3.cd(2);
339  hLRtotSoverSplusB -> Draw();
340  c3.Print(fname,"Landscape");
341 
342  TCanvas c4("c4","",1);
343  hEffvsPur -> Draw("AL*");
344  c4.Print(fname,"Landscape");
345 
346  TCanvas c5("c5","",1);
347  hLRValvsPur -> Draw("AL*");
348  c5.Print(fname,"Landscape");
349 
350  TCanvas c6("c6","",1);
351  hLRValvsEff -> Draw("AL*");
352  c6.Print(fname,"Landscape");
353  }
354  c.Print(fname + "]","landscape");
355 }
std::vector< TH1F * > hObsSoverSplusB
std::vector< TH1F * > hObsS
std::vector< TH1F * > hObsB
TGraph * hLRValvsEff
std::vector< TH2F * > hObsCorr
TGraph * hLRValvsPur
string fname
main script
std::vector< TF1 * > fObsSoverSplusB
TH1F * hLRtotSoverSplusB
void LRHelpFunctions::storeToROOTfile ( const TString &  fname)

Definition at line 256 of file LRHelpFunctions.cc.

References constructPurity, fLRtotSoverSplusB, fObsSoverSplusB, hEffvsPur, hLRtotB, hLRtotS, hLRtotSoverSplusB, hLRValvsEff, hLRValvsPur, hObsB, hObsCorr, hObsS, hObsSoverSplusB, and python.connectstrParser::o.

256  {
257  TFile fOut(fname,"RECREATE");
258  fOut.cd();
259  for(size_t o=0; o<hObsS.size(); o++){
260  hObsS[o] -> Write();
261  hObsB[o] -> Write();
262  hObsSoverSplusB[o] -> Write();
263  fObsSoverSplusB[o] -> Write();
264  }
265  int hIndex = 0;
266  for(size_t o=0; o<hObsS.size(); o++) {
267  for(size_t o2=o+1; o2<hObsS.size(); o2++) {
268  hObsCorr[hIndex] -> Write();
269  ++ hIndex;
270  }
271  }
272  if(constructPurity){
273  hLRtotS -> Write();
274  hLRtotB -> Write();
275  hLRtotSoverSplusB -> Write();
276  fLRtotSoverSplusB -> Write();
277  hEffvsPur -> Write();
278  hLRValvsPur -> Write();
279  hLRValvsEff -> Write();
280  }
281  fOut.cd();
282  fOut.Write();
283  fOut.Close();
284 }
std::vector< TH1F * > hObsSoverSplusB
std::vector< TH1F * > hObsS
std::vector< TH1F * > hObsB
TGraph * hLRValvsEff
std::vector< TH2F * > hObsCorr
TGraph * hLRValvsPur
string fname
main script
std::vector< TF1 * > fObsSoverSplusB
TH1F * hLRtotSoverSplusB

Member Data Documentation

bool LRHelpFunctions::constructPurity
private
TF1* LRHelpFunctions::fLRtotSoverSplusB
private
std::vector<TF1*> LRHelpFunctions::fObsSoverSplusB
private
TGraph* LRHelpFunctions::hEffvsPur
private
TH1F * LRHelpFunctions::hLRtotB
private
TH1F* LRHelpFunctions::hLRtotS
private
TH1F * LRHelpFunctions::hLRtotSoverSplusB
private
TGraph * LRHelpFunctions::hLRValvsEff
private
TGraph * LRHelpFunctions::hLRValvsPur
private
std::vector<TH1F*> LRHelpFunctions::hObsB
private
std::vector<TH2F*> LRHelpFunctions::hObsCorr
private
std::vector<TH1F*> LRHelpFunctions::hObsS
private
std::vector<TH1F*> LRHelpFunctions::hObsSoverSplusB
private
std::vector<int> LRHelpFunctions::obsNumbers
private

Definition at line 79 of file LRHelpFunctions.h.

Referenced by calcPtdrLRval(), LRHelpFunctions(), and readObsHistsAndFits().