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 Member Functions | Private Attributes
EffPurFromHistos Class Reference

#include <EffPurFromHistos.h>

Public Member Functions

void compute ()
 
FlavourHistograms< double > * discriminatorCutEfficScan () const
 
FlavourHistograms< double > * discriminatorNoCutEffic () const
 
 EffPurFromHistos (const std::string &ext, TH1F *h_d, TH1F *h_u, TH1F *h_s, TH1F *h_c, TH1F *h_b, TH1F *h_g, TH1F *h_ni, TH1F *h_dus, TH1F *h_dusg, const std::string &label, const bool &mc, int nBin=100, double startO=0.005, double endO=1.005)
 
 EffPurFromHistos (const FlavourHistograms< double > *dDiscriminatorFC, const std::string &label, const bool &mc, int nBin=100, double startO=0.005, double endO=1.005)
 
void epsPlot (const std::string &name)
 
TH1F * getEffFlavVsBEff_b ()
 
TH1F * getEffFlavVsBEff_c ()
 
TH1F * getEffFlavVsBEff_d ()
 
TH1F * getEffFlavVsBEff_dus ()
 
TH1F * getEffFlavVsBEff_dusg ()
 
TH1F * getEffFlavVsBEff_g ()
 
TH1F * getEffFlavVsBEff_ni ()
 
TH1F * getEffFlavVsBEff_s ()
 
TH1F * getEffFlavVsBEff_u ()
 
void plot (TPad *theCanvas=0)
 
void plot (const std::string &name, const std::string &ext)
 
void psPlot (const std::string &name)
 
 ~EffPurFromHistos ()
 

Private Member Functions

void check ()
 

Private Attributes

FlavourHistograms< double > * discrCutEfficScan
 
FlavourHistograms< double > * discrNoCutEffic
 
MonitorElementEffFlavVsBEff_b
 
MonitorElementEffFlavVsBEff_c
 
MonitorElementEffFlavVsBEff_d
 
MonitorElementEffFlavVsBEff_dus
 
MonitorElementEffFlavVsBEff_dusg
 
MonitorElementEffFlavVsBEff_g
 
MonitorElementEffFlavVsBEff_ni
 
MonitorElementEffFlavVsBEff_s
 
MonitorElementEffFlavVsBEff_u
 
TH1F * effVersusDiscr_b
 
TH1F * effVersusDiscr_c
 
TH1F * effVersusDiscr_d
 
TH1F * effVersusDiscr_dus
 
TH1F * effVersusDiscr_dusg
 
TH1F * effVersusDiscr_g
 
TH1F * effVersusDiscr_ni
 
TH1F * effVersusDiscr_s
 
TH1F * effVersusDiscr_u
 
double endOutput
 
bool fromDiscriminatorDistr
 
std::string histoExtension
 
std::string label_
 
bool mcPlots_
 
int nBinOutput
 
double startOutput
 

Detailed Description

Definition at line 14 of file EffPurFromHistos.h.

Constructor & Destructor Documentation

EffPurFromHistos::EffPurFromHistos ( const std::string &  ext,
TH1F *  h_d,
TH1F *  h_u,
TH1F *  h_s,
TH1F *  h_c,
TH1F *  h_b,
TH1F *  h_g,
TH1F *  h_ni,
TH1F *  h_dus,
TH1F *  h_dusg,
const std::string &  label,
const bool &  mc,
int  nBin = 100,
double  startO = 0.005,
double  endO = 1.005 
)

Definition at line 20 of file EffPurFromHistos.cc.

References check().

22  :
23  //BTagPlotPrintC(),
28  effVersusDiscr_dusg(h_dusg), nBinOutput(nBin), startOutput(startO),
29  endOutput(endO), mcPlots_(mc), label_(label)
30 {
31  // consistency check
32  check();
33 }
std::string histoExtension
EffPurFromHistos::EffPurFromHistos ( const FlavourHistograms< double > *  dDiscriminatorFC,
const std::string &  label,
const bool &  mc,
int  nBin = 100,
double  startO = 0.005,
double  endO = 1.005 
)

Definition at line 36 of file EffPurFromHistos.cc.

References FlavourHistograms< T >::baseNameDescription(), FlavourHistograms< T >::baseNameTitle(), FlavourHistograms< T >::divide(), FlavourHistograms< T >::getHistoVector(), FlavourHistograms< T >::histo_b(), FlavourHistograms< T >::histo_c(), FlavourHistograms< T >::histo_d(), FlavourHistograms< T >::histo_dus(), FlavourHistograms< T >::histo_dusg(), FlavourHistograms< T >::histo_g(), FlavourHistograms< T >::histo_ni(), FlavourHistograms< T >::histo_s(), FlavourHistograms< T >::histo_u(), diffTwoXMLs::label, FlavourHistograms< T >::lowerBound(), FlavourHistograms< T >::nBins(), FlavourHistograms< T >::SetMinimum(), mathSSE::sqrt(), and FlavourHistograms< T >::upperBound().

37  :
38  fromDiscriminatorDistr(true), nBinOutput(nBin), startOutput(startO), endOutput(endO), mcPlots_(mc), label_(label){
39  histoExtension = "_"+dDiscriminatorFC->baseNameTitle();
40 
41 
43  "totalEntries" + histoExtension, "Total Entries: " + dDiscriminatorFC->baseNameDescription(),
44  dDiscriminatorFC->nBins(), dDiscriminatorFC->lowerBound(),
45  dDiscriminatorFC->upperBound(), false, true, false, "b", false, label, mcPlots_ );
46 
47  // conditional discriminator cut for efficiency histos
48 
50  "effVsDiscrCut" + histoExtension, "Eff. vs Disc. Cut: " + dDiscriminatorFC->baseNameDescription(),
51  dDiscriminatorFC->nBins(), dDiscriminatorFC->lowerBound(),
52  dDiscriminatorFC->upperBound(), false, true, false, "b", false, label , mcPlots_ );
53  discrCutEfficScan->SetMinimum(1E-4);
54  if (mcPlots_){
55 
56  effVersusDiscr_d = discrCutEfficScan->histo_d ();
57  effVersusDiscr_u = discrCutEfficScan->histo_u ();
58  effVersusDiscr_s = discrCutEfficScan->histo_s ();
59  effVersusDiscr_c = discrCutEfficScan->histo_c ();
60  effVersusDiscr_b = discrCutEfficScan->histo_b ();
61  effVersusDiscr_g = discrCutEfficScan->histo_g ();
62  effVersusDiscr_ni = discrCutEfficScan->histo_ni ();
63  effVersusDiscr_dus = discrCutEfficScan->histo_dus ();
64  effVersusDiscr_dusg = discrCutEfficScan->histo_dusg();
65 
66 
67 
68  effVersusDiscr_d->SetXTitle ( "Discriminant" );
69  effVersusDiscr_d->GetXaxis()->SetTitleOffset ( 0.75 );
70  effVersusDiscr_u->SetXTitle ( "Discriminant" );
71  effVersusDiscr_u->GetXaxis()->SetTitleOffset ( 0.75 );
72  effVersusDiscr_s->SetXTitle ( "Discriminant" );
73  effVersusDiscr_s->GetXaxis()->SetTitleOffset ( 0.75 );
74  effVersusDiscr_c->SetXTitle ( "Discriminant" );
75  effVersusDiscr_c->GetXaxis()->SetTitleOffset ( 0.75 );
76  effVersusDiscr_b->SetXTitle ( "Discriminant" );
77  effVersusDiscr_b->GetXaxis()->SetTitleOffset ( 0.75 );
78  effVersusDiscr_g->SetXTitle ( "Discriminant" );
79  effVersusDiscr_g->GetXaxis()->SetTitleOffset ( 0.75 );
80  effVersusDiscr_ni->SetXTitle ( "Discriminant" );
81  effVersusDiscr_ni->GetXaxis()->SetTitleOffset ( 0.75 );
82  effVersusDiscr_dus->SetXTitle ( "Discriminant" );
83  effVersusDiscr_dus->GetXaxis()->SetTitleOffset ( 0.75 );
84  effVersusDiscr_dusg->SetXTitle ( "Discriminant" );
85  effVersusDiscr_dusg->GetXaxis()->SetTitleOffset ( 0.75 );
86  }else{
87  effVersusDiscr_d = 0;
88  effVersusDiscr_u = 0;
89  effVersusDiscr_s = 0;
90  effVersusDiscr_c = 0;
91  effVersusDiscr_b = 0;
92  effVersusDiscr_g = 0;
96 }
97 
98  // discr. for computation
99  vector<TH1F*> discrCfHistos = dDiscriminatorFC->getHistoVector();
100 
101  // discr no cut
102  vector<TH1F*> discrNoCutHistos = discrNoCutEffic->getHistoVector();
103 
104  // discr no cut
105  vector<TH1F*> discrCutHistos = discrCutEfficScan->getHistoVector();
106 
107  const int& dimHistos = discrCfHistos.size(); // they all have the same size
108 
109  // DISCR-CUT LOOP:
110  // fill the histos for eff-pur computations by scanning the discriminatorFC histogram
111 
112  // better to loop over bins -> discrCut no longer needed
113  const int& nBins = dDiscriminatorFC->nBins();
114 
115  // loop over flavours
116  for ( int iFlav = 0; iFlav < dimHistos; iFlav++ ) {
117  if (discrCfHistos[iFlav] == 0) continue;
118  discrNoCutHistos[iFlav]->SetXTitle ( "Discriminant" );
119  discrNoCutHistos[iFlav]->GetXaxis()->SetTitleOffset ( 0.75 );
120 
121  // In Root histos, bin counting starts at 1 to nBins.
122  // bin 0 is the underflow, and nBins+1 is the overflow.
123  const double& nJetsFlav = discrCfHistos[iFlav]->GetEntries ();
124  double sum = discrCfHistos[iFlav]->GetBinContent( nBins+1 ); //+1 to get the overflow.
125 
126  for ( int iDiscr = nBins; iDiscr > 0 ; --iDiscr ) {
127  // fill all jets into NoCut histo
128  discrNoCutHistos[iFlav]->SetBinContent ( iDiscr, nJetsFlav );
129  discrNoCutHistos[iFlav]->SetBinError ( iDiscr, sqrt(nJetsFlav) );
130  sum += discrCfHistos[iFlav]->GetBinContent( iDiscr );
131  discrCutHistos[iFlav]->SetBinContent ( iDiscr, sum );
132  discrCutHistos[iFlav]->SetBinError ( iDiscr, sqrt(sum) );
133  }
134  }
135 
136 
137  // divide to get efficiency vs. discriminator cut from absolute numbers
138  discrCutEfficScan->divide ( *discrNoCutEffic ); // does: histos including discriminator cut / flat histo
139 }
double lowerBound() const
void SetMinimum(const double &min)
TH1F * histo_c() const
FlavourHistograms< double > * discrNoCutEffic
std::string baseNameDescription() const
TH1F * histo_dus() const
std::vector< TH1F * > getHistoVector() const
TH1F * histo_b() const
TH1F * histo_g() const
T sqrt(T t)
Definition: SSEVec.h:46
std::string baseNameTitle() const
TH1F * histo_d() const
FlavourHistograms< double > * discrCutEfficScan
std::string histoExtension
TH1F * histo_s() const
TH1F * histo_ni() const
void divide(const FlavourHistograms< T > &bHD) const
TH1F * histo_u() const
double upperBound() const
TH1F * histo_dusg() const
EffPurFromHistos::~EffPurFromHistos ( )

Definition at line 142 of file EffPurFromHistos.cc.

142  {
143  /* delete EffFlavVsBEff_d ;
144  delete EffFlavVsBEff_u ;
145  delete EffFlavVsBEff_s ;
146  delete EffFlavVsBEff_c ;
147  delete EffFlavVsBEff_b ;
148  delete EffFlavVsBEff_g ;
149  delete EffFlavVsBEff_ni ;
150  delete EffFlavVsBEff_dus ;
151  delete EffFlavVsBEff_dusg;
152  if ( fromDiscriminatorDistr) {
153  delete discrNoCutEffic;
154  delete discrCutEfficScan;
155  }*/
156 }

Member Function Documentation

void EffPurFromHistos::check ( )
private

Definition at line 313 of file EffPurFromHistos.cc.

References effVersusDiscr_b, effVersusDiscr_c, effVersusDiscr_d, effVersusDiscr_dus, effVersusDiscr_dusg, effVersusDiscr_g, effVersusDiscr_ni, effVersusDiscr_s, effVersusDiscr_u, and edm::hlt::Exception.

Referenced by EffPurFromHistos().

313  {
314  // number of bins
315  const int& nBins_d = effVersusDiscr_d -> GetNbinsX();
316  const int& nBins_u = effVersusDiscr_u -> GetNbinsX();
317  const int& nBins_s = effVersusDiscr_s -> GetNbinsX();
318  const int& nBins_c = effVersusDiscr_c -> GetNbinsX();
319  const int& nBins_b = effVersusDiscr_b -> GetNbinsX();
320  const int& nBins_g = effVersusDiscr_g -> GetNbinsX();
321  const int& nBins_ni = effVersusDiscr_ni -> GetNbinsX();
322  const int& nBins_dus = effVersusDiscr_dus -> GetNbinsX();
323  const int& nBins_dusg = effVersusDiscr_dusg -> GetNbinsX();
324 
325  const bool& lNBins =
326  ( nBins_d == nBins_u &&
327  nBins_d == nBins_s &&
328  nBins_d == nBins_c &&
329  nBins_d == nBins_b &&
330  nBins_d == nBins_g &&
331  nBins_d == nBins_ni &&
332  nBins_d == nBins_dus &&
333  nBins_d == nBins_dusg );
334 
335  if ( !lNBins ) {
336  throw cms::Exception("Configuration")
337  << "Input histograms do not all have the same number of bins!\n";
338  }
339 
340 
341  // start
342  const float& sBin_d = effVersusDiscr_d -> GetBinCenter(1);
343  const float& sBin_u = effVersusDiscr_u -> GetBinCenter(1);
344  const float& sBin_s = effVersusDiscr_s -> GetBinCenter(1);
345  const float& sBin_c = effVersusDiscr_c -> GetBinCenter(1);
346  const float& sBin_b = effVersusDiscr_b -> GetBinCenter(1);
347  const float& sBin_g = effVersusDiscr_g -> GetBinCenter(1);
348  const float& sBin_ni = effVersusDiscr_ni -> GetBinCenter(1);
349  const float& sBin_dus = effVersusDiscr_dus -> GetBinCenter(1);
350  const float& sBin_dusg = effVersusDiscr_dusg -> GetBinCenter(1);
351 
352  const bool& lSBin =
353  ( sBin_d == sBin_u &&
354  sBin_d == sBin_s &&
355  sBin_d == sBin_c &&
356  sBin_d == sBin_b &&
357  sBin_d == sBin_g &&
358  sBin_d == sBin_ni &&
359  sBin_d == sBin_dus &&
360  sBin_d == sBin_dusg );
361 
362  if ( !lSBin ) {
363  throw cms::Exception("Configuration")
364  << "EffPurFromHistos::check() : Input histograms do not all have the same start bin!\n";
365  }
366 
367 
368  // end
369  const float& eBin_d = effVersusDiscr_d -> GetBinCenter( nBins_d - 1 );
370  const float& eBin_u = effVersusDiscr_u -> GetBinCenter( nBins_d - 1 );
371  const float& eBin_s = effVersusDiscr_s -> GetBinCenter( nBins_d - 1 );
372  const float& eBin_c = effVersusDiscr_c -> GetBinCenter( nBins_d - 1 );
373  const float& eBin_b = effVersusDiscr_b -> GetBinCenter( nBins_d - 1 );
374  const float& eBin_g = effVersusDiscr_g -> GetBinCenter( nBins_d - 1 );
375  const float& eBin_ni = effVersusDiscr_ni -> GetBinCenter( nBins_d - 1 );
376  const float& eBin_dus = effVersusDiscr_dus -> GetBinCenter( nBins_d - 1 );
377  const float& eBin_dusg = effVersusDiscr_dusg -> GetBinCenter( nBins_d - 1 );
378 
379  const bool& lEBin =
380  ( eBin_d == eBin_u &&
381  eBin_d == eBin_s &&
382  eBin_d == eBin_c &&
383  eBin_d == eBin_b &&
384  eBin_d == eBin_g &&
385  eBin_d == eBin_ni &&
386  eBin_d == eBin_dus &&
387  eBin_d == eBin_dusg );
388 
389  if ( !lEBin ) {
390  throw cms::Exception("Configuration")
391  << "EffPurFromHistos::check() : Input histograms do not all have the same end bin!\n";
392  }
393 }
void EffPurFromHistos::compute ( )

Definition at line 396 of file EffPurFromHistos.cc.

References HistoProviderDQM::book1D(), EffFlavVsBEff_b, EffFlavVsBEff_c, EffFlavVsBEff_d, EffFlavVsBEff_dus, EffFlavVsBEff_dusg, EffFlavVsBEff_g, EffFlavVsBEff_ni, EffFlavVsBEff_s, EffFlavVsBEff_u, effVersusDiscr_b, effVersusDiscr_c, effVersusDiscr_d, effVersusDiscr_dus, effVersusDiscr_dusg, effVersusDiscr_g, effVersusDiscr_ni, effVersusDiscr_s, effVersusDiscr_u, endOutput, HcalObjRepresent::Fill(), RecoBTag::findBinClosestYValue(), MonitorElement::getTH1F(), histoExtension, label_, mcPlots_, nBinOutput, and startOutput.

Referenced by JetTagPlotter::finalize().

397 {
398  if (!mcPlots_) {
399 
400  EffFlavVsBEff_d = 0;
401  EffFlavVsBEff_u = 0;
402  EffFlavVsBEff_s = 0;
403  EffFlavVsBEff_c = 0;
404  EffFlavVsBEff_b = 0;
405  EffFlavVsBEff_g = 0;
406  EffFlavVsBEff_ni = 0;
407  EffFlavVsBEff_dus = 0;
408  EffFlavVsBEff_dusg = 0;
409 
410  return;
411 
412  }
413 
414 
415  // to have shorter names ......
416  const std::string & hE = histoExtension;
417  const std::string & hB = "FlavEffVsBEff_";
418 
419 
420  // create histograms from base name and extension as given from user
421  // BINNING MUST BE IDENTICAL FOR ALL OF THEM!!
422  HistoProviderDQM prov("Btag",label_);
423  EffFlavVsBEff_d = (prov.book1D ( hB + "D" + hE , hB + "D" + hE , nBinOutput , startOutput , endOutput ));
424  EffFlavVsBEff_u = (prov.book1D ( hB + "U" + hE , hB + "U" + hE , nBinOutput , startOutput , endOutput )) ;
425  EffFlavVsBEff_s = (prov.book1D ( hB + "S" + hE , hB + "S" + hE , nBinOutput , startOutput , endOutput )) ;
426  EffFlavVsBEff_c = (prov.book1D ( hB + "C" + hE , hB + "C" + hE , nBinOutput , startOutput , endOutput )) ;
427  EffFlavVsBEff_b = (prov.book1D ( hB + "B" + hE , hB + "B" + hE , nBinOutput , startOutput , endOutput )) ;
428  EffFlavVsBEff_g = (prov.book1D ( hB + "G" + hE , hB + "G" + hE , nBinOutput , startOutput , endOutput )) ;
429  EffFlavVsBEff_ni = (prov.book1D ( hB + "NI" + hE , hB + "NI" + hE , nBinOutput , startOutput , endOutput )) ;
430  EffFlavVsBEff_dus = (prov.book1D ( hB + "DUS" + hE , hB + "DUS" + hE , nBinOutput , startOutput , endOutput )) ;
431  EffFlavVsBEff_dusg = (prov.book1D ( hB + "DUSG" + hE , hB + "DUSG" + hE , nBinOutput , startOutput , endOutput )) ;
432 
433  EffFlavVsBEff_d->getTH1F()->SetXTitle ( "b-jet efficiency" );
434  EffFlavVsBEff_d->getTH1F()->SetYTitle ( "non b-jet efficiency" );
435  EffFlavVsBEff_d->getTH1F()->GetXaxis()->SetTitleOffset ( 0.75 );
436  EffFlavVsBEff_d->getTH1F()->GetYaxis()->SetTitleOffset ( 0.75 );
437  EffFlavVsBEff_u->getTH1F()->SetXTitle ( "b-jet efficiency" );
438  EffFlavVsBEff_u->getTH1F()->SetYTitle ( "non b-jet efficiency" );
439  EffFlavVsBEff_u->getTH1F()->GetXaxis()->SetTitleOffset ( 0.75 );
440  EffFlavVsBEff_u->getTH1F()->GetYaxis()->SetTitleOffset ( 0.75 );
441  EffFlavVsBEff_s->getTH1F()->SetXTitle ( "b-jet efficiency" );
442  EffFlavVsBEff_s->getTH1F()->SetYTitle ( "non b-jet efficiency" );
443  EffFlavVsBEff_s->getTH1F()->GetXaxis()->SetTitleOffset ( 0.75 );
444  EffFlavVsBEff_s->getTH1F()->GetYaxis()->SetTitleOffset ( 0.75 );
445  EffFlavVsBEff_c->getTH1F()->SetXTitle ( "b-jet efficiency" );
446  EffFlavVsBEff_c->getTH1F()->SetYTitle ( "non b-jet efficiency" );
447  EffFlavVsBEff_s->getTH1F()->GetXaxis()->SetTitleOffset ( 0.75 );
448  EffFlavVsBEff_s->getTH1F()->GetYaxis()->SetTitleOffset ( 0.75 );
449  EffFlavVsBEff_b->getTH1F()->SetXTitle ( "b-jet efficiency" );
450  EffFlavVsBEff_b->getTH1F()->SetYTitle ( "b-jet efficiency" );
451  EffFlavVsBEff_b->getTH1F()->GetXaxis()->SetTitleOffset ( 0.75 );
452  EffFlavVsBEff_b->getTH1F()->GetYaxis()->SetTitleOffset ( 0.75 );
453  EffFlavVsBEff_g->getTH1F()->SetXTitle ( "b-jet efficiency" );
454  EffFlavVsBEff_g->getTH1F()->SetYTitle ( "non b-jet efficiency" );
455  EffFlavVsBEff_g->getTH1F()->GetXaxis()->SetTitleOffset ( 0.75 );
456  EffFlavVsBEff_g->getTH1F()->GetYaxis()->SetTitleOffset ( 0.75 );
457  EffFlavVsBEff_ni->getTH1F()->SetXTitle ( "b-jet efficiency" );
458  EffFlavVsBEff_ni->getTH1F()->SetYTitle ( "non b-jet efficiency" );
459  EffFlavVsBEff_ni->getTH1F()->GetXaxis()->SetTitleOffset ( 0.75 );
460  EffFlavVsBEff_ni->getTH1F()->GetYaxis()->SetTitleOffset ( 0.75 );
461  EffFlavVsBEff_dus->getTH1F()->SetXTitle ( "b-jet efficiency" );
462  EffFlavVsBEff_dus->getTH1F()->SetYTitle ( "non b-jet efficiency" );
463  EffFlavVsBEff_dus->getTH1F()->GetXaxis()->SetTitleOffset ( 0.75 );
464  EffFlavVsBEff_dus->getTH1F()->GetYaxis()->SetTitleOffset ( 0.75 );
465  EffFlavVsBEff_dusg->getTH1F()->SetXTitle ( "b-jet efficiency" );
466  EffFlavVsBEff_dusg->getTH1F()->SetYTitle ( "non b-jet efficiency" );
467  EffFlavVsBEff_dusg->getTH1F()->GetXaxis()->SetTitleOffset ( 0.75 );
468  EffFlavVsBEff_dusg->getTH1F()->GetYaxis()->SetTitleOffset ( 0.75 );
469 
470 
471  // loop over eff. vs. discriminator cut b-histo and look in which bin the closest entry is;
472  // use fact that eff decreases monotonously
473 
474  // any of the histos to be created can be taken here:
475  MonitorElement * EffFlavVsBEff = EffFlavVsBEff_b;
476 
477  const int& nBinB = EffFlavVsBEff->getTH1F()->GetNbinsX();
478 
479  for ( int iBinB = 1; iBinB <= nBinB; iBinB++ ) { // loop over the bins on the x-axis of the histograms to be filled
480 
481  const float& effBBinWidth = EffFlavVsBEff->getTH1F()->GetBinWidth ( iBinB );
482  const float& effBMid = EffFlavVsBEff->getTH1F()->GetBinCenter ( iBinB ); // middle of b-efficiency bin
483  const float& effBLeft = effBMid - 0.5*effBBinWidth; // left edge of bin
484  const float& effBRight = effBMid + 0.5*effBBinWidth; // right edge of bin
485  // find the corresponding bin in the efficiency versus discriminator cut histo: closest one in efficiency
486  const int& binClosest = findBinClosestYValue ( effVersusDiscr_b , effBMid , effBLeft , effBRight );
487  const bool& binFound = ( binClosest > 0 ) ;
488  //
489  if ( binFound ) {
490  // fill the histos
491  EffFlavVsBEff_d -> Fill ( effBMid , effVersusDiscr_d ->GetBinContent ( binClosest ) );
492  EffFlavVsBEff_u -> Fill ( effBMid , effVersusDiscr_u ->GetBinContent ( binClosest ) );
493  EffFlavVsBEff_s -> Fill ( effBMid , effVersusDiscr_s ->GetBinContent ( binClosest ) );
494  EffFlavVsBEff_c -> Fill ( effBMid , effVersusDiscr_c ->GetBinContent ( binClosest ) );
495  EffFlavVsBEff_b -> Fill ( effBMid , effVersusDiscr_b ->GetBinContent ( binClosest ) );
496  EffFlavVsBEff_g -> Fill ( effBMid , effVersusDiscr_g ->GetBinContent ( binClosest ) );
497  EffFlavVsBEff_ni -> Fill ( effBMid , effVersusDiscr_ni ->GetBinContent ( binClosest ) );
498  EffFlavVsBEff_dus -> Fill ( effBMid , effVersusDiscr_dus ->GetBinContent ( binClosest ) );
499  EffFlavVsBEff_dusg -> Fill ( effBMid , effVersusDiscr_dusg->GetBinContent ( binClosest ) );
500 
501  EffFlavVsBEff_d ->getTH1F() -> SetBinError ( iBinB , effVersusDiscr_d ->GetBinError ( binClosest ) );
502  EffFlavVsBEff_u ->getTH1F() -> SetBinError ( iBinB , effVersusDiscr_u ->GetBinError ( binClosest ) );
503  EffFlavVsBEff_s ->getTH1F() -> SetBinError ( iBinB , effVersusDiscr_s ->GetBinError ( binClosest ) );
504  EffFlavVsBEff_c ->getTH1F() -> SetBinError ( iBinB , effVersusDiscr_c ->GetBinError ( binClosest ) );
505  EffFlavVsBEff_b ->getTH1F() -> SetBinError ( iBinB , effVersusDiscr_b ->GetBinError ( binClosest ) );
506  EffFlavVsBEff_g ->getTH1F() -> SetBinError ( iBinB , effVersusDiscr_g ->GetBinError ( binClosest ) );
507  EffFlavVsBEff_ni ->getTH1F() -> SetBinError ( iBinB , effVersusDiscr_ni ->GetBinError ( binClosest ) );
508  EffFlavVsBEff_dus->getTH1F() -> SetBinError ( iBinB , effVersusDiscr_dus ->GetBinError ( binClosest ) );
509  EffFlavVsBEff_dusg->getTH1F() -> SetBinError ( iBinB , effVersusDiscr_dusg->GetBinError ( binClosest ) );
510  }
511  else {
512  //CW cout << "Did not find right bin for b-efficiency : " << effBMid << endl;
513  }
514 
515  }
516 
517 }
MonitorElement * EffFlavVsBEff_g
MonitorElement * EffFlavVsBEff_dusg
MonitorElement * EffFlavVsBEff_u
int findBinClosestYValue(const TH1F *, const float &yVal, const float &yLow, const float &yHigh)
Definition: Tools.cc:233
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
MonitorElement * EffFlavVsBEff_c
MonitorElement * EffFlavVsBEff_s
std::string histoExtension
TH1F * getTH1F(void) const
MonitorElement * EffFlavVsBEff_d
MonitorElement * EffFlavVsBEff_b
MonitorElement * EffFlavVsBEff_ni
MonitorElement * EffFlavVsBEff_dus
FlavourHistograms<double>* EffPurFromHistos::discriminatorCutEfficScan ( ) const
inline
FlavourHistograms<double>* EffPurFromHistos::discriminatorNoCutEffic ( ) const
inline
void EffPurFromHistos::epsPlot ( const std::string &  name)

Definition at line 161 of file EffPurFromHistos.cc.

References discrCutEfficScan, discrNoCutEffic, FlavourHistograms< T >::epsPlot(), fromDiscriminatorDistr, and plot().

Referenced by JetTagPlotter::epsPlot().

162 {
163  if ( fromDiscriminatorDistr) {
166  }
167  plot(name, ".eps");
168 }
FlavourHistograms< double > * discrNoCutEffic
void epsPlot(const std::string &name)
void plot(TPad *theCanvas=0)
FlavourHistograms< double > * discrCutEfficScan
TH1F* EffPurFromHistos::getEffFlavVsBEff_b ( )
inline

Definition at line 40 of file EffPurFromHistos.h.

References EffFlavVsBEff_b, and MonitorElement::getTH1F().

Referenced by BTagDifferentialPlot::fillHisto().

40 { return EffFlavVsBEff_b ->getTH1F() ; };
TH1F * getTH1F(void) const
MonitorElement * EffFlavVsBEff_b
TH1F* EffPurFromHistos::getEffFlavVsBEff_c ( )
inline

Definition at line 39 of file EffPurFromHistos.h.

References EffFlavVsBEff_c, and MonitorElement::getTH1F().

Referenced by BTagDifferentialPlot::fillHisto().

39 { return EffFlavVsBEff_c ->getTH1F() ; };
MonitorElement * EffFlavVsBEff_c
TH1F * getTH1F(void) const
TH1F* EffPurFromHistos::getEffFlavVsBEff_d ( )
inline

Definition at line 35 of file EffPurFromHistos.h.

References EffFlavVsBEff_d, and MonitorElement::getTH1F().

Referenced by BTagDifferentialPlot::fillHisto().

35  {
36  return EffFlavVsBEff_d->getTH1F() ; };
TH1F * getTH1F(void) const
MonitorElement * EffFlavVsBEff_d
TH1F* EffPurFromHistos::getEffFlavVsBEff_dus ( )
inline

Definition at line 43 of file EffPurFromHistos.h.

References EffFlavVsBEff_dus, and MonitorElement::getTH1F().

Referenced by BTagDifferentialPlot::fillHisto().

43 { return EffFlavVsBEff_dus ->getTH1F() ; };
TH1F * getTH1F(void) const
MonitorElement * EffFlavVsBEff_dus
TH1F* EffPurFromHistos::getEffFlavVsBEff_dusg ( )
inline

Definition at line 44 of file EffPurFromHistos.h.

References EffFlavVsBEff_dusg, and MonitorElement::getTH1F().

Referenced by BTagDifferentialPlot::fillHisto().

44 { return EffFlavVsBEff_dusg ->getTH1F(); };
MonitorElement * EffFlavVsBEff_dusg
TH1F * getTH1F(void) const
TH1F* EffPurFromHistos::getEffFlavVsBEff_g ( )
inline

Definition at line 41 of file EffPurFromHistos.h.

References EffFlavVsBEff_g, and MonitorElement::getTH1F().

Referenced by BTagDifferentialPlot::fillHisto().

41 { return EffFlavVsBEff_g ->getTH1F() ; };
MonitorElement * EffFlavVsBEff_g
TH1F * getTH1F(void) const
TH1F* EffPurFromHistos::getEffFlavVsBEff_ni ( )
inline

Definition at line 42 of file EffPurFromHistos.h.

References EffFlavVsBEff_ni, and MonitorElement::getTH1F().

Referenced by BTagDifferentialPlot::fillHisto().

42 { return EffFlavVsBEff_ni ->getTH1F() ; };
TH1F * getTH1F(void) const
MonitorElement * EffFlavVsBEff_ni
TH1F* EffPurFromHistos::getEffFlavVsBEff_s ( )
inline

Definition at line 38 of file EffPurFromHistos.h.

References EffFlavVsBEff_s, and MonitorElement::getTH1F().

Referenced by BTagDifferentialPlot::fillHisto().

38 { return EffFlavVsBEff_s ->getTH1F() ; };
MonitorElement * EffFlavVsBEff_s
TH1F * getTH1F(void) const
TH1F* EffPurFromHistos::getEffFlavVsBEff_u ( )
inline

Definition at line 37 of file EffPurFromHistos.h.

References EffFlavVsBEff_u, and MonitorElement::getTH1F().

Referenced by BTagDifferentialPlot::fillHisto().

37 { return EffFlavVsBEff_u->getTH1F() ; };
MonitorElement * EffFlavVsBEff_u
TH1F * getTH1F(void) const
void EffPurFromHistos::plot ( TPad *  theCanvas = 0)

Definition at line 183 of file EffPurFromHistos.cc.

References alignCSCRings::e, EffFlavVsBEff_b, EffFlavVsBEff_c, EffFlavVsBEff_d, EffFlavVsBEff_dus, EffFlavVsBEff_dusg, EffFlavVsBEff_g, EffFlavVsBEff_ni, EffFlavVsBEff_s, EffFlavVsBEff_u, MonitorElement::getTH1F(), and plotscripts::setTDRStyle().

Referenced by epsPlot(), plot(), TrackCountingTagPlotter::psPlot(), TrackIPTagPlotter::psPlot(), TrackProbabilityTagPlotter::psPlot(), JetTagPlotter::psPlot(), and psPlot().

183  {
184 
185 //fixme:
186  bool btppNI = false;
187  bool btppColour = true;
188 
189 // if ( !btppTitle ) gStyle->SetOptTitle ( 0 );
190  setTDRStyle()->cd();
191 
192  if (plotCanvas)
193  plotCanvas->cd();
194 
195  gPad->UseCurrentStyle();
196  gPad->SetFillColor ( 0 );
197  gPad->SetLogy ( 1 );
198  gPad->SetGridx ( 1 );
199  gPad->SetGridy ( 1 );
200 
201  int col_c ;
202  int col_g ;
203  int col_dus;
204  int col_ni ;
205 
206  int mStyle_c ;
207  int mStyle_g ;
208  int mStyle_dus;
209  int mStyle_ni ;
210 
211  // marker size (same for all)
212  float mSize = gPad->GetWh() * gPad->GetHNDC() / 500.; //1.2;
213 
214  if ( btppColour ) {
215  col_c = 6;
216  col_g = 3; // g in green
217  col_dus = 4; // uds in blue
218  col_ni = 5; // ni in ??
219  mStyle_c = 20;
220  mStyle_g = 20;
221  mStyle_dus = 20;
222  mStyle_ni = 20;
223  }
224  else {
225  col_c = 1;
226  col_g = 1;
227  col_dus = 1;
228  col_ni = 1;
229  mStyle_c = 22;
230  mStyle_g = 29;
231  mStyle_dus = 20;
232  mStyle_ni = 27;
233  }
234 
235 
236  // for the moment: plot c,dus,g
237  EffFlavVsBEff_dus ->getTH1F()->GetXaxis()->SetTitle ( "b-jet efficiency" );
238  EffFlavVsBEff_dus ->getTH1F()->GetYaxis()->SetTitle ( "non b-jet efficiency" );
239  EffFlavVsBEff_dus ->getTH1F()->GetYaxis()->SetTitleOffset ( 0.25 );
240  EffFlavVsBEff_dus ->getTH1F()->SetMaximum ( 1.1 );
241  EffFlavVsBEff_dus ->getTH1F()->SetMinimum ( 1.e-5 );
242  EffFlavVsBEff_dus ->getTH1F()->SetMarkerColor ( col_dus );
243  EffFlavVsBEff_dus ->getTH1F()->SetLineColor ( col_dus );
244  EffFlavVsBEff_dus ->getTH1F()->SetMarkerSize ( mSize );
245  EffFlavVsBEff_dus ->getTH1F()->SetMarkerStyle ( mStyle_dus );
246  EffFlavVsBEff_dus ->getTH1F()->SetStats ( false );
247  EffFlavVsBEff_dus ->getTH1F()->Draw("pe");
248 
249  EffFlavVsBEff_g ->getTH1F()->SetMarkerColor ( col_g );
250  EffFlavVsBEff_g ->getTH1F()->SetLineColor ( col_g );
251  EffFlavVsBEff_g ->getTH1F()->SetMarkerSize ( mSize );
252  EffFlavVsBEff_g ->getTH1F()->SetMarkerStyle ( mStyle_g );
253  EffFlavVsBEff_g ->getTH1F()->SetStats ( false );
254  EffFlavVsBEff_g ->getTH1F()->Draw("peSame");
255 
256  EffFlavVsBEff_c ->getTH1F()->SetMarkerColor ( col_c );
257  EffFlavVsBEff_c ->getTH1F()->SetLineColor ( col_c );
258  EffFlavVsBEff_c ->getTH1F()->SetMarkerSize ( mSize );
259  EffFlavVsBEff_c ->getTH1F()->SetMarkerStyle ( mStyle_c );
260  EffFlavVsBEff_c ->getTH1F()->SetStats ( false );
261  EffFlavVsBEff_c ->getTH1F()->Draw("peSame");
262 
263  EffFlavVsBEff_d ->getTH1F()-> SetMinimum(0.01);
264  EffFlavVsBEff_u ->getTH1F()-> SetMinimum(0.01);
265  EffFlavVsBEff_s ->getTH1F()-> SetMinimum(0.01);
266  EffFlavVsBEff_c ->getTH1F()-> SetMinimum(0.01);
267  EffFlavVsBEff_b ->getTH1F()-> SetMinimum(0.01);
268  EffFlavVsBEff_g ->getTH1F()-> SetMinimum(0.01);
269  EffFlavVsBEff_ni ->getTH1F()-> SetMinimum(0.01);
270  EffFlavVsBEff_dus ->getTH1F()-> SetMinimum(0.01);
271  EffFlavVsBEff_dusg ->getTH1F()-> SetMinimum(0.01);
272 
273  // plot separately u,d and s
274 // EffFlavVsBEff_d ->GetXaxis()->SetTitle ( "b-jet efficiency" );
275 // EffFlavVsBEff_d ->GetYaxis()->SetTitle ( "non b-jet efficiency" );
276 // EffFlavVsBEff_d ->GetYaxis()->SetTitleOffset ( 1.25 );
277 // EffFlavVsBEff_d ->SetMaximum ( 1.1 );
278 // EffFlavVsBEff_d ->SetMinimum ( 1.e-5 );
279 // EffFlavVsBEff_d ->SetMarkerColor ( col_dus );
280 // EffFlavVsBEff_d ->SetLineColor ( col_dus );
281 // EffFlavVsBEff_d ->SetMarkerSize ( mSize );
282 // EffFlavVsBEff_d ->SetMarkerStyle ( mStyle_dus );
283 // EffFlavVsBEff_d ->SetStats ( false );
284 // EffFlavVsBEff_d ->Draw("pe");
285 //
286 // EffFlavVsBEff_u ->SetMarkerColor ( col_g );
287 // EffFlavVsBEff_u ->SetLineColor ( col_g );
288 // EffFlavVsBEff_u ->SetMarkerSize ( mSize );
289 // EffFlavVsBEff_u ->SetMarkerStyle ( mStyle_g );
290 // EffFlavVsBEff_u ->SetStats ( false );
291 // EffFlavVsBEff_u ->Draw("peSame");
292 //
293 // EffFlavVsBEff_s ->SetMarkerColor ( col_c );
294 // EffFlavVsBEff_s ->SetLineColor ( col_c );
295 // EffFlavVsBEff_s ->SetMarkerSize ( mSize );
296 // EffFlavVsBEff_s ->SetMarkerStyle ( mStyle_c );
297 // EffFlavVsBEff_s ->SetStats ( false );
298 // EffFlavVsBEff_s ->Draw("peSame");
299 
300  // only if asked: NI
301  if ( btppNI ) {
302  EffFlavVsBEff_ni ->getTH1F()->SetMarkerColor ( col_ni );
303  EffFlavVsBEff_ni ->getTH1F()->SetLineColor ( col_ni );
304  EffFlavVsBEff_ni ->getTH1F()->SetMarkerSize ( mSize );
305  EffFlavVsBEff_ni ->getTH1F()->SetMarkerStyle ( mStyle_ni );
306  EffFlavVsBEff_ni ->getTH1F()->SetStats ( false );
307  EffFlavVsBEff_ni ->getTH1F()->Draw("peSame");
308  }
309 
310 }
MonitorElement * EffFlavVsBEff_g
MonitorElement * EffFlavVsBEff_dusg
def setTDRStyle
Definition: plotscripts.py:87
MonitorElement * EffFlavVsBEff_u
MonitorElement * EffFlavVsBEff_c
MonitorElement * EffFlavVsBEff_s
TH1F * getTH1F(void) const
MonitorElement * EffFlavVsBEff_d
MonitorElement * EffFlavVsBEff_b
MonitorElement * EffFlavVsBEff_ni
MonitorElement * EffFlavVsBEff_dus
void EffPurFromHistos::plot ( const std::string &  name,
const std::string &  ext 
)

Definition at line 175 of file EffPurFromHistos.cc.

References histoExtension, and plot().

176 {
177  TCanvas tc (("FlavEffVsBEff" +histoExtension).c_str() ,
178  ("Flavour misidentification vs. b-tagging efficiency " + histoExtension).c_str());
179  plot(&tc);
180  tc.Print((name + "FlavEffVsBEff" + histoExtension + ext).c_str());
181 }
void plot(TPad *theCanvas=0)
std::string histoExtension
void EffPurFromHistos::psPlot ( const std::string &  name)

Definition at line 170 of file EffPurFromHistos.cc.

References plot().

171 {
172  plot(name, ".ps");
173 }
void plot(TPad *theCanvas=0)

Member Data Documentation

FlavourHistograms<double> * EffPurFromHistos::discrCutEfficScan
private

Definition at line 72 of file EffPurFromHistos.h.

Referenced by discriminatorCutEfficScan(), and epsPlot().

FlavourHistograms<double>* EffPurFromHistos::discrNoCutEffic
private

Definition at line 72 of file EffPurFromHistos.h.

Referenced by discriminatorNoCutEffic(), and epsPlot().

MonitorElement* EffPurFromHistos::EffFlavVsBEff_b
private

Definition at line 103 of file EffPurFromHistos.h.

Referenced by compute(), getEffFlavVsBEff_b(), and plot().

MonitorElement* EffPurFromHistos::EffFlavVsBEff_c
private

Definition at line 102 of file EffPurFromHistos.h.

Referenced by compute(), getEffFlavVsBEff_c(), and plot().

MonitorElement* EffPurFromHistos::EffFlavVsBEff_d
private

Definition at line 99 of file EffPurFromHistos.h.

Referenced by compute(), getEffFlavVsBEff_d(), and plot().

MonitorElement* EffPurFromHistos::EffFlavVsBEff_dus
private

Definition at line 106 of file EffPurFromHistos.h.

Referenced by compute(), getEffFlavVsBEff_dus(), and plot().

MonitorElement* EffPurFromHistos::EffFlavVsBEff_dusg
private

Definition at line 107 of file EffPurFromHistos.h.

Referenced by compute(), getEffFlavVsBEff_dusg(), and plot().

MonitorElement* EffPurFromHistos::EffFlavVsBEff_g
private

Definition at line 104 of file EffPurFromHistos.h.

Referenced by compute(), getEffFlavVsBEff_g(), and plot().

MonitorElement* EffPurFromHistos::EffFlavVsBEff_ni
private

Definition at line 105 of file EffPurFromHistos.h.

Referenced by compute(), getEffFlavVsBEff_ni(), and plot().

MonitorElement* EffPurFromHistos::EffFlavVsBEff_s
private

Definition at line 101 of file EffPurFromHistos.h.

Referenced by compute(), getEffFlavVsBEff_s(), and plot().

MonitorElement* EffPurFromHistos::EffFlavVsBEff_u
private

Definition at line 100 of file EffPurFromHistos.h.

Referenced by compute(), getEffFlavVsBEff_u(), and plot().

TH1F* EffPurFromHistos::effVersusDiscr_b
private

Definition at line 82 of file EffPurFromHistos.h.

Referenced by check(), and compute().

TH1F* EffPurFromHistos::effVersusDiscr_c
private

Definition at line 81 of file EffPurFromHistos.h.

Referenced by check(), and compute().

TH1F* EffPurFromHistos::effVersusDiscr_d
private

Definition at line 78 of file EffPurFromHistos.h.

Referenced by check(), and compute().

TH1F* EffPurFromHistos::effVersusDiscr_dus
private

Definition at line 85 of file EffPurFromHistos.h.

Referenced by check(), and compute().

TH1F* EffPurFromHistos::effVersusDiscr_dusg
private

Definition at line 86 of file EffPurFromHistos.h.

Referenced by check(), and compute().

TH1F* EffPurFromHistos::effVersusDiscr_g
private

Definition at line 83 of file EffPurFromHistos.h.

Referenced by check(), and compute().

TH1F* EffPurFromHistos::effVersusDiscr_ni
private

Definition at line 84 of file EffPurFromHistos.h.

Referenced by check(), and compute().

TH1F* EffPurFromHistos::effVersusDiscr_s
private

Definition at line 80 of file EffPurFromHistos.h.

Referenced by check(), and compute().

TH1F* EffPurFromHistos::effVersusDiscr_u
private

Definition at line 79 of file EffPurFromHistos.h.

Referenced by check(), and compute().

double EffPurFromHistos::endOutput
private

Definition at line 94 of file EffPurFromHistos.h.

Referenced by compute().

bool EffPurFromHistos::fromDiscriminatorDistr
private

Definition at line 66 of file EffPurFromHistos.h.

Referenced by epsPlot().

std::string EffPurFromHistos::histoExtension
private

Definition at line 70 of file EffPurFromHistos.h.

Referenced by compute(), and plot().

std::string EffPurFromHistos::label_
private
bool EffPurFromHistos::mcPlots_
private

Definition at line 96 of file EffPurFromHistos.h.

Referenced by compute().

int EffPurFromHistos::nBinOutput
private

Definition at line 92 of file EffPurFromHistos.h.

Referenced by compute().

double EffPurFromHistos::startOutput
private

Definition at line 93 of file EffPurFromHistos.h.

Referenced by compute().