CMS 3D CMS Logo

List of all members | Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes
FlavourHistograms< T > Class Template Reference

#include <FlavourHistorgrams.h>

Inheritance diagram for FlavourHistograms< T >:
TrackIPHistograms< T >

Public Member Functions

int * arrayDimension () const
 
std::string baseNameDescription () const
 
std::string baseNameTitle () const
 
void divide (const FlavourHistograms< T > &bHD)
 
void epsPlot (const std::string &name)
 
void fill (const int &flavour, const T &variable) const
 
void fill (const int &flavour, const T &variable, const T &w) const
 
void fill (const int &flavour, const T *variable) const
 
 FlavourHistograms (const std::string &baseNameTitle_, const std::string &baseNameDescription_, const int &nBins_, const double &lowerBound_, const double &upperBound_, const std::string &plotFirst_, const std::string &folder, const unsigned int &mc, DQMStore::IGetter &iget)
 
 FlavourHistograms (const std::string &baseNameTitle_, const std::string &baseNameDescription_, const int &nBins_, const double &lowerBound_, const double &upperBound_, const bool &statistics_, const bool &plotLog_, const bool &plotNormalized_, const std::string &plotFirst_, const std::string &folder, const unsigned int &mc, DQMStore::IBooker &ibook)
 
std::vector< TH1F * > getHistoVector () const
 
TH1F * histo_all () const
 
TH1F * histo_b () const
 
TH1F * histo_c () const
 
TH1F * histo_d () const
 
TH1F * histo_dus () const
 
TH1F * histo_dusg () const
 
TH1F * histo_g () const
 
TH1F * histo_ni () const
 
TH1F * histo_pu () const
 
TH1F * histo_s () const
 
TH1F * histo_u () const
 
int indexToPlot () const
 
double lowerBound () const
 
int maxDimension () const
 
int nBins () const
 
void plot (TPad *theCanvas=0)
 
std::string plotFirst () const
 
bool plotLog () const
 
bool plotNormalized () const
 
void setEfficiencyFlag ()
 
void SetMaximum (const double &max)
 
void SetMinimum (const double &min)
 
void settitle (const char *title)
 
bool statistics () const
 
double upperBound () const
 
virtual ~FlavourHistograms ()
 

Protected Member Functions

double ClopperPearsonUnc (double num, double den)
 
void ComputeEfficiency (TH1F *num, TH1F *den, int bin)
 
void fillVariable (const int &flavour, const T &var, const T &w) const
 

Protected Attributes

int * theArrayDimension
 
std::string theBaseNameDescription
 
std::string theBaseNameTitle
 
TCanvas * theCanvas
 
MonitorElementtheHisto_all
 
MonitorElementtheHisto_b
 
MonitorElementtheHisto_c
 
MonitorElementtheHisto_d
 
MonitorElementtheHisto_dus
 
MonitorElementtheHisto_dusg
 
MonitorElementtheHisto_g
 
MonitorElementtheHisto_ni
 
MonitorElementtheHisto_pu
 
MonitorElementtheHisto_s
 
MonitorElementtheHisto_u
 
int theIndexToPlot
 
double theLowerBound
 
double theMax
 
int theMaxDimension
 
double theMin
 
int theNBins
 
std::string thePlotFirst
 
bool thePlotLog
 
bool thePlotNormalized
 
bool theStatistics
 
double theUpperBound
 

Private Member Functions

 FlavourHistograms ()
 

Private Attributes

unsigned int mcPlots_
 

Detailed Description

template<class T>
class FlavourHistograms< T >

Definition at line 29 of file FlavourHistorgrams.h.

Constructor & Destructor Documentation

template<class T >
FlavourHistograms< T >::FlavourHistograms ( const std::string &  baseNameTitle_,
const std::string &  baseNameDescription_,
const int &  nBins_,
const double &  lowerBound_,
const double &  upperBound_,
const std::string &  plotFirst_,
const std::string &  folder,
const unsigned int &  mc,
DQMStore::IGetter iget 
)

Definition at line 151 of file FlavourHistorgrams.h.

References gather_cfg::cout, DQMStore::IGetter::get(), FlavourHistograms< T >::mcPlots_, nullptr, FlavourHistograms< T >::theArrayDimension, FlavourHistograms< T >::theBaseNameTitle, FlavourHistograms< T >::theCanvas, FlavourHistograms< T >::theHisto_all, FlavourHistograms< T >::theHisto_b, FlavourHistograms< T >::theHisto_c, FlavourHistograms< T >::theHisto_d, FlavourHistograms< T >::theHisto_dus, FlavourHistograms< T >::theHisto_dusg, FlavourHistograms< T >::theHisto_g, FlavourHistograms< T >::theHisto_ni, FlavourHistograms< T >::theHisto_pu, FlavourHistograms< T >::theHisto_s, FlavourHistograms< T >::theHisto_u, and FlavourHistograms< T >::thePlotFirst.

154  :
155  theMaxDimension(-1), theIndexToPlot(-1), theBaseNameTitle ( baseNameTitle_ ) , theBaseNameDescription ( baseNameDescription_ ) ,
156  theNBins ( nBins_ ) , theLowerBound ( lowerBound_ ) , theUpperBound ( upperBound_ ) ,
157  theStatistics ( false ) , thePlotLog ( false ) , thePlotNormalized ( false ) ,
158  thePlotFirst ( plotFirst_ ), theMin(-1.), theMax(-1.), mcPlots_(mc)
159 {
160  // defaults for array dimensions
162 
163  // check plot order string
164  if ( thePlotFirst == "l" || thePlotFirst == "c" || thePlotFirst == "b" ) {
165  // OK
166  }
167  else {
168  // not correct: print warning and set default (l)
169  std::cout << "FlavourHistograms::FlavourHistograms : thePlotFirst was not correct : " << thePlotFirst << std::endl ;
170  std::cout << "FlavourHistograms::FlavourHistograms : Set it to default value (l)! " << std::endl ;
171  thePlotFirst = "l" ;
172  }
173 
174  if(mcPlots_%2==0) theHisto_all = iget.get("Btag/" + folder + "/" +theBaseNameTitle + "ALL" ) ;
175  else theHisto_all = nullptr;
176  if (mcPlots_) {
177  if (mcPlots_>2 ) {
178  theHisto_d = iget.get("Btag/" + folder + "/" +theBaseNameTitle + "D" ) ;
179  theHisto_u = iget.get("Btag/" + folder + "/" +theBaseNameTitle + "U" ) ;
180  theHisto_s = iget.get("Btag/" + folder + "/" +theBaseNameTitle + "S" ) ;
181  theHisto_g = iget.get("Btag/" + folder + "/" +theBaseNameTitle + "G" ) ;
182  theHisto_dus = iget.get("Btag/" + folder + "/" +theBaseNameTitle + "DUS" ) ;
183  }
184  else{
185  theHisto_d = nullptr;
186  theHisto_u = nullptr;
187  theHisto_s = nullptr;
188  theHisto_g = nullptr;
189  theHisto_dus = nullptr;
190  }
191  theHisto_c = iget.get("Btag/" + folder + "/" +theBaseNameTitle + "C" ) ;
192  theHisto_b = iget.get("Btag/" + folder + "/" +theBaseNameTitle + "B" ) ;
193  theHisto_ni = iget.get("Btag/" + folder + "/" +theBaseNameTitle + "NI" ) ;
194  theHisto_dusg = iget.get("Btag/" + folder + "/" +theBaseNameTitle + "DUSG") ;
195  theHisto_pu = iget.get("Btag/" + folder + "/" +theBaseNameTitle + "PU") ;
196  }
197  else{
198  theHisto_d = nullptr;
199  theHisto_u = nullptr;
200  theHisto_s = nullptr;
201  theHisto_c = nullptr;
202  theHisto_b = nullptr;
203  theHisto_g = nullptr;
204  theHisto_ni = nullptr;
205  theHisto_dus = nullptr;
206  theHisto_dusg = nullptr;
207  theHisto_pu = nullptr;
208  }
209  // defaults for other data members
210  theCanvas = nullptr ;
211 }
MonitorElement * theHisto_g
MonitorElement * theHisto_c
MonitorElement * get(const std::string &path)
Definition: DQMStore.cc:307
MonitorElement * theHisto_u
MonitorElement * theHisto_all
#define nullptr
MonitorElement * theHisto_pu
MonitorElement * theHisto_d
MonitorElement * theHisto_dusg
std::string theBaseNameDescription
MonitorElement * theHisto_b
std::string theBaseNameTitle
MonitorElement * theHisto_s
MonitorElement * theHisto_dus
MonitorElement * theHisto_ni
template<class T >
FlavourHistograms< T >::FlavourHistograms ( const std::string &  baseNameTitle_,
const std::string &  baseNameDescription_,
const int &  nBins_,
const double &  lowerBound_,
const double &  upperBound_,
const bool &  statistics_,
const bool &  plotLog_,
const bool &  plotNormalized_,
const std::string &  plotFirst_,
const std::string &  folder,
const unsigned int &  mc,
DQMStore::IBooker ibook 
)

Definition at line 214 of file FlavourHistorgrams.h.

References HistoProviderDQM::book1D(), gather_cfg::cout, MonitorElement::getTH1F(), FlavourHistograms< T >::mcPlots_, nullptr, FlavourHistograms< T >::theArrayDimension, FlavourHistograms< T >::theBaseNameDescription, FlavourHistograms< T >::theBaseNameTitle, FlavourHistograms< T >::theCanvas, FlavourHistograms< T >::theHisto_all, FlavourHistograms< T >::theHisto_b, FlavourHistograms< T >::theHisto_c, FlavourHistograms< T >::theHisto_d, FlavourHistograms< T >::theHisto_dus, FlavourHistograms< T >::theHisto_dusg, FlavourHistograms< T >::theHisto_g, FlavourHistograms< T >::theHisto_ni, FlavourHistograms< T >::theHisto_pu, FlavourHistograms< T >::theHisto_s, FlavourHistograms< T >::theHisto_u, FlavourHistograms< T >::theLowerBound, FlavourHistograms< T >::theNBins, FlavourHistograms< T >::thePlotFirst, FlavourHistograms< T >::theStatistics, and FlavourHistograms< T >::theUpperBound.

218  :
219  theMaxDimension(-1), theIndexToPlot(-1), theBaseNameTitle ( baseNameTitle_ ) , theBaseNameDescription ( baseNameDescription_ ) ,
220  theNBins ( nBins_ ) , theLowerBound ( lowerBound_ ) , theUpperBound ( upperBound_ ) ,
221  theStatistics ( statistics_ ) , thePlotLog ( plotLog_ ) , thePlotNormalized ( plotNormalized_ ) ,
222  thePlotFirst ( plotFirst_ ), theMin(-1.), theMax(-1.), mcPlots_(mc)
223 {
224  // defaults for array dimensions
226 
227  // check plot order string
228  if ( thePlotFirst == "l" || thePlotFirst == "c" || thePlotFirst == "b" ) {
229  // OK
230  }
231  else {
232  // not correct: print warning and set default (l)
233  std::cout << "FlavourHistograms::FlavourHistograms : thePlotFirst was not correct : " << thePlotFirst << std::endl ;
234  std::cout << "FlavourHistograms::FlavourHistograms : Set it to default value (l)! " << std::endl ;
235  thePlotFirst = "l" ;
236  }
237 
238  // book histos
239  HistoProviderDQM prov("Btag",folder,ibook);
240  if(mcPlots_%2==0) theHisto_all = (prov.book1D( theBaseNameTitle + "ALL" , theBaseNameDescription + " all jets" , theNBins , theLowerBound , theUpperBound )) ;
241  else theHisto_all = nullptr;
242  if (mcPlots_) {
243  if(mcPlots_>2) {
244  theHisto_d = (prov.book1D ( theBaseNameTitle + "D" , theBaseNameDescription + " d-jets" , theNBins , theLowerBound , theUpperBound )) ;
245  theHisto_u = (prov.book1D ( theBaseNameTitle + "U" , theBaseNameDescription + " u-jets" , theNBins , theLowerBound , theUpperBound )) ;
246  theHisto_s = (prov.book1D ( theBaseNameTitle + "S" , theBaseNameDescription + " s-jets" , theNBins , theLowerBound , theUpperBound )) ;
247  theHisto_g = (prov.book1D ( theBaseNameTitle + "G" , theBaseNameDescription + " g-jets" , theNBins , theLowerBound , theUpperBound )) ;
248  theHisto_dus = (prov.book1D ( theBaseNameTitle + "DUS" , theBaseNameDescription + " dus-jets" , theNBins , theLowerBound , theUpperBound )) ;
249  }
250  else{
251  theHisto_d = nullptr;
252  theHisto_u = nullptr;
253  theHisto_s = nullptr;
254  theHisto_g = nullptr;
255  theHisto_dus = nullptr;
256  }
257  theHisto_c = (prov.book1D ( theBaseNameTitle + "C" , theBaseNameDescription + " c-jets" , theNBins , theLowerBound , theUpperBound )) ;
258  theHisto_b = (prov.book1D ( theBaseNameTitle + "B" , theBaseNameDescription + " b-jets" , theNBins , theLowerBound , theUpperBound )) ;
259  theHisto_ni = (prov.book1D ( theBaseNameTitle + "NI" , theBaseNameDescription + " ni-jets" , theNBins , theLowerBound , theUpperBound )) ;
260  theHisto_dusg = (prov.book1D ( theBaseNameTitle + "DUSG" , theBaseNameDescription + " dusg-jets" , theNBins , theLowerBound , theUpperBound )) ;
261  theHisto_pu = (prov.book1D ( theBaseNameTitle + "PU" , theBaseNameDescription + " pu-jets" , theNBins , theLowerBound , theUpperBound )) ;
262  }else{
263  theHisto_d = nullptr;
264  theHisto_u = nullptr;
265  theHisto_s = nullptr;
266  theHisto_c = nullptr;
267  theHisto_b = nullptr;
268  theHisto_g = nullptr;
269  theHisto_ni = nullptr;
270  theHisto_dus = nullptr;
271  theHisto_dusg = nullptr;
272  theHisto_pu = nullptr;
273  }
274 
275  // statistics if requested
276  if ( theStatistics ) {
277  if(theHisto_all) theHisto_all ->getTH1F()->Sumw2() ;
278  if (mcPlots_ ) {
279  if (mcPlots_>2 ) {
280  theHisto_d ->getTH1F()->Sumw2() ;
281  theHisto_u ->getTH1F()->Sumw2() ;
282  theHisto_s ->getTH1F()->Sumw2() ;
283  theHisto_g ->getTH1F()->Sumw2() ;
284  theHisto_dus ->getTH1F()->Sumw2() ;
285  }
286  theHisto_c ->getTH1F()->Sumw2() ;
287  theHisto_b ->getTH1F()->Sumw2() ;
288  theHisto_ni ->getTH1F()->Sumw2() ;
289  theHisto_dusg->getTH1F()->Sumw2() ;
290  theHisto_pu ->getTH1F()->Sumw2() ;
291  }
292  }
293  // defaults for other data members
294  theCanvas = nullptr ;
295 }
MonitorElement * theHisto_g
MonitorElement * theHisto_c
TH1F * getTH1F() const
MonitorElement * theHisto_u
MonitorElement * theHisto_all
#define nullptr
MonitorElement * theHisto_pu
MonitorElement * theHisto_d
MonitorElement * theHisto_dusg
std::string theBaseNameDescription
MonitorElement * theHisto_b
std::string theBaseNameTitle
MonitorElement * theHisto_s
MonitorElement * theHisto_dus
MonitorElement * theHisto_ni
template<class T >
FlavourHistograms< T >::~FlavourHistograms ( )
virtual

Definition at line 299 of file FlavourHistorgrams.h.

References FlavourHistograms< T >::theCanvas.

299  {
300  // delete the canvas*/
301  delete theCanvas ;
302 }
template<class T>
FlavourHistograms< T >::FlavourHistograms ( )
inlineprivate

Definition at line 144 of file FlavourHistorgrams.h.

144 {}

Member Function Documentation

template<class T>
int* FlavourHistograms< T >::arrayDimension ( ) const
inline
template<class T>
std::string FlavourHistograms< T >::baseNameDescription ( ) const
inline
template<class T>
std::string FlavourHistograms< T >::baseNameTitle ( ) const
inline

Definition at line 68 of file FlavourHistorgrams.h.

References FlavourHistograms< T >::theBaseNameTitle.

Referenced by EffPurFromHistos::EffPurFromHistos().

68 { return theBaseNameTitle ; }
std::string theBaseNameTitle
template<class T >
double FlavourHistograms< T >::ClopperPearsonUnc ( double  num,
double  den 
)
protected

Definition at line 511 of file FlavourHistorgrams.h.

References SiStripPI::max.

Referenced by FlavourHistograms< T >::ComputeEfficiency(), and FlavourHistograms< T >::histo_pu().

511  {
512  double effVal = num / den;
513  double errLo = TEfficiency::ClopperPearson(static_cast<int>(den), static_cast<int>(num), 0.683, false);
514  double errUp = TEfficiency::ClopperPearson(static_cast<int>(den), static_cast<int>(num), 0.683, true);
515  return std::max(effVal - errLo, errUp - effVal);
516 }
template<class T >
void FlavourHistograms< T >::ComputeEfficiency ( TH1F *  num,
TH1F *  den,
int  bin 
)
protected

Definition at line 519 of file FlavourHistorgrams.h.

References FlavourHistograms< T >::ClopperPearsonUnc().

Referenced by FlavourHistograms< T >::divide(), and FlavourHistograms< T >::histo_pu().

519  {
520  double effVal = 1.;
521  double errVal = 0.;
522  double numVal = num->GetBinContent(bin);
523  double denVal = den->GetBinContent(bin);
524  if (denVal > 0 && numVal <= denVal) {
525  effVal = numVal / denVal;
526  errVal = ClopperPearsonUnc(numVal, denVal);
527  }
528  num->SetBinContent(bin, effVal);
529  num->SetBinError(bin, errVal);
530 }
double ClopperPearsonUnc(double num, double den)
bin
set the eta bin as selection string.
template<class T >
void FlavourHistograms< T >::divide ( const FlavourHistograms< T > &  bHD)

Definition at line 533 of file FlavourHistorgrams.h.

References stringResolutionProvider_cfi::bin, FlavourHistograms< T >::ComputeEfficiency(), FlavourHistograms< T >::histo_all(), 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_pu(), FlavourHistograms< T >::histo_s(), FlavourHistograms< T >::histo_u(), FlavourHistograms< T >::mcPlots_, FlavourHistograms< T >::theHisto_all, FlavourHistograms< T >::theHisto_b, FlavourHistograms< T >::theHisto_c, FlavourHistograms< T >::theHisto_d, FlavourHistograms< T >::theHisto_dus, FlavourHistograms< T >::theHisto_dusg, FlavourHistograms< T >::theHisto_g, FlavourHistograms< T >::theHisto_ni, FlavourHistograms< T >::theHisto_pu, FlavourHistograms< T >::theHisto_s, FlavourHistograms< T >::theHisto_u, and FlavourHistograms< T >::theNBins.

Referenced by python.rootplot.utilities.Hist::__div__().

533  {
534  for(int bin = 0; bin < theNBins+2; bin++){
535  if(theHisto_all) ComputeEfficiency(theHisto_all ->getTH1F(), bHD.histo_all(), bin);
536  if (mcPlots_) {
537  if (mcPlots_>2 ) {
538  ComputeEfficiency(theHisto_d ->getTH1F(), bHD.histo_d (), bin) ;
539  ComputeEfficiency(theHisto_u ->getTH1F(), bHD.histo_u (), bin) ;
540  ComputeEfficiency(theHisto_s ->getTH1F(), bHD.histo_s (), bin) ;
541  ComputeEfficiency(theHisto_g ->getTH1F(), bHD.histo_g (), bin) ;
542  ComputeEfficiency(theHisto_dus ->getTH1F(), bHD.histo_dus (), bin) ;
543  }
544  ComputeEfficiency(theHisto_c ->getTH1F(), bHD.histo_c (), bin) ;
545  ComputeEfficiency(theHisto_b ->getTH1F(), bHD.histo_b (), bin) ;
546  ComputeEfficiency(theHisto_ni ->getTH1F(), bHD.histo_ni (), bin) ;
547  ComputeEfficiency(theHisto_dusg ->getTH1F(), bHD.histo_dusg(), bin) ;
548  ComputeEfficiency(theHisto_pu ->getTH1F(), bHD.histo_pu (), bin) ;
549  }
550  }
551 }
MonitorElement * theHisto_g
MonitorElement * theHisto_c
TH1F * histo_all() const
MonitorElement * theHisto_u
TH1F * histo_c() const
MonitorElement * theHisto_all
TH1F * histo_dus() const
void ComputeEfficiency(TH1F *num, TH1F *den, int bin)
TH1F * histo_b() const
TH1F * histo_g() const
MonitorElement * theHisto_pu
TH1F * histo_d() const
MonitorElement * theHisto_d
MonitorElement * theHisto_dusg
bin
set the eta bin as selection string.
MonitorElement * theHisto_b
TH1F * histo_s() const
TH1F * histo_ni() const
TH1F * histo_u() const
MonitorElement * theHisto_s
TH1F * histo_dusg() const
MonitorElement * theHisto_dus
TH1F * histo_pu() const
MonitorElement * theHisto_ni
template<class T >
void FlavourHistograms< T >::epsPlot ( const std::string &  name)

Definition at line 502 of file FlavourHistorgrams.h.

References FlavourHistograms< T >::plot(), FlavourHistograms< T >::theBaseNameDescription, and FlavourHistograms< T >::theBaseNameTitle.

503 {
504  TCanvas tc(theBaseNameTitle.c_str() , theBaseNameDescription.c_str());
505 
506  plot(&tc);
507  tc.Print((name + theBaseNameTitle + ".eps").c_str());
508 }
void plot(TPad *theCanvas=0)
std::string theBaseNameDescription
std::string theBaseNameTitle
template<class T >
void FlavourHistograms< T >::fill ( const int &  flavour,
const T variable 
) const

Definition at line 306 of file FlavourHistorgrams.h.

References FlavourHistograms< T >::fillVariable().

Referenced by TrackIPHistograms< T >::fill().

307 {
308  // For single variables and arrays (for arrays only a single index can be filled)
309  fillVariable ( flavour , variable, 1.) ;
310 }
void fillVariable(const int &flavour, const T &var, const T &w) const
template<class T >
void FlavourHistograms< T >::fill ( const int &  flavour,
const T variable,
const T w 
) const

Definition at line 313 of file FlavourHistorgrams.h.

References FlavourHistograms< T >::fillVariable().

314 {
315  // For single variables and arrays (for arrays only a single index can be filled)
317 }
const double w
Definition: UKUtility.cc:23
void fillVariable(const int &flavour, const T &var, const T &w) const
template<class T >
void FlavourHistograms< T >::fill ( const int &  flavour,
const T variable 
) const

Definition at line 320 of file FlavourHistorgrams.h.

References FlavourHistograms< T >::fillVariable(), mps_fire::i, FlavourHistograms< T >::theArrayDimension, FlavourHistograms< T >::theIndexToPlot, and FlavourHistograms< T >::theMaxDimension.

321 {
322  if ( theArrayDimension == 0 ) {
323  // single variable
324  fillVariable ( flavour , *variable, 1.) ;
325  } else {
326  // array
328  //
329  for ( int i = 0 ; i != iMax ; ++i ) {
330  // check if only one index to be plotted (<0: switched off -> plot all)
331  if ( ( theIndexToPlot < 0 ) || ( i == theIndexToPlot ) ) {
332  fillVariable ( flavour , *(variable+i) , 1.) ;
333  }
334  }
335 
336  // if single index to be filled but not enough entries: fill 0.0 (convention!)
337  if ( theIndexToPlot >= iMax ) {
338  // cout << "==>> The index to be filled is too big -> fill 0.0 : " << theBaseNameTitle << " : " << theIndexToPlot << " >= " << iMax << endl ;
339  const T& theZero = static_cast<T> ( 0.0 ) ;
340  fillVariable ( flavour , theZero , 1.) ;
341  }
342  }
343 }
void fillVariable(const int &flavour, const T &var, const T &w) const
long double T
template<class T >
void FlavourHistograms< T >::fillVariable ( const int &  flavour,
const T var,
const T w 
) const
protected

Definition at line 574 of file FlavourHistorgrams.h.

References MonitorElement::Fill(), FlavourHistograms< T >::mcPlots_, FlavourHistograms< T >::theBaseNameDescription, FlavourHistograms< T >::theHisto_all, FlavourHistograms< T >::theHisto_b, FlavourHistograms< T >::theHisto_c, FlavourHistograms< T >::theHisto_d, FlavourHistograms< T >::theHisto_dus, FlavourHistograms< T >::theHisto_dusg, FlavourHistograms< T >::theHisto_g, FlavourHistograms< T >::theHisto_ni, FlavourHistograms< T >::theHisto_pu, FlavourHistograms< T >::theHisto_s, and FlavourHistograms< T >::theHisto_u.

Referenced by FlavourHistograms< T >::fill(), and FlavourHistograms< T >::histo_pu().

574  {
575 
576  // all, except for the Jet Multiplicity which is not filled for each jets but for each events
577  if( (theBaseNameDescription != "Jet Multiplicity" || flavour == -1) && theHisto_all) theHisto_all->Fill ( var ,w) ;
578 
579  // flavour specific
580  if (!mcPlots_ || (theBaseNameDescription == "Jet Multiplicity" && flavour == -1)) return;
581 
582  switch(flavour) {
583  case 1:
584  if (mcPlots_>2 ) {
585  theHisto_d->Fill( var ,w);
586  if(theBaseNameDescription != "Jet Multiplicity") theHisto_dus->Fill( var ,w);
587  }
588  if(theBaseNameDescription != "Jet Multiplicity") theHisto_dusg->Fill( var ,w);
589  return;
590  case 2:
591  if (mcPlots_>2 ) {
592  theHisto_u->Fill( var ,w);
593  if(theBaseNameDescription != "Jet Multiplicity") theHisto_dus->Fill( var ,w);
594  }
595  if(theBaseNameDescription != "Jet Multiplicity") theHisto_dusg->Fill( var ,w);
596  return;
597  case 3:
598  if (mcPlots_>2 ) {
599  theHisto_s->Fill( var ,w);
600  if(theBaseNameDescription != "Jet Multiplicity") theHisto_dus->Fill( var ,w);
601  }
602  if(theBaseNameDescription != "Jet Multiplicity") theHisto_dusg->Fill( var ,w);
603  return;
604  case 4:
605  theHisto_c->Fill( var ,w);
606  return;
607  case 5:
608  theHisto_b->Fill( var ,w);
609  return;
610  case 21:
611  if (mcPlots_>2 ) theHisto_g->Fill( var ,w);
612  if(theBaseNameDescription != "Jet Multiplicity") theHisto_dusg->Fill( var ,w);
613  return;
614  case 123:
615  if (mcPlots_>2 && theBaseNameDescription == "Jet Multiplicity") theHisto_dus->Fill( var ,w);
616  return;
617  case 12321:
618  if (theBaseNameDescription == "Jet Multiplicity") theHisto_dusg->Fill( var ,w);
619  return;
620  case 20:
621  theHisto_pu->Fill( var ,w);
622  return;
623  default:
624  theHisto_ni->Fill( var ,w);
625  return;
626  }
627 }
MonitorElement * theHisto_g
MonitorElement * theHisto_c
const double w
Definition: UKUtility.cc:23
MonitorElement * theHisto_u
MonitorElement * theHisto_all
void Fill(long long x)
MonitorElement * theHisto_pu
MonitorElement * theHisto_d
MonitorElement * theHisto_dusg
std::string theBaseNameDescription
MonitorElement * theHisto_b
MonitorElement * theHisto_s
MonitorElement * theHisto_dus
MonitorElement * theHisto_ni
template<class T >
std::vector< TH1F * > FlavourHistograms< T >::getHistoVector ( ) const

Definition at line 630 of file FlavourHistorgrams.h.

References MonitorElement::getTH1F(), FlavourHistograms< T >::mcPlots_, FlavourHistograms< T >::theHisto_all, FlavourHistograms< T >::theHisto_b, FlavourHistograms< T >::theHisto_c, FlavourHistograms< T >::theHisto_d, FlavourHistograms< T >::theHisto_dus, FlavourHistograms< T >::theHisto_dusg, FlavourHistograms< T >::theHisto_g, FlavourHistograms< T >::theHisto_ni, FlavourHistograms< T >::theHisto_pu, FlavourHistograms< T >::theHisto_s, and FlavourHistograms< T >::theHisto_u.

Referenced by EffPurFromHistos::EffPurFromHistos(), and FlavourHistograms< T >::histo_pu().

631 {
632  std::vector<TH1F*> histoVector;
633  if(theHisto_all) histoVector.push_back ( theHisto_all->getTH1F() );
634  if (mcPlots_) {
635  if (mcPlots_>2 ) {
636  histoVector.push_back ( theHisto_d->getTH1F() );
637  histoVector.push_back ( theHisto_u->getTH1F() );
638  histoVector.push_back ( theHisto_s->getTH1F() );
639  histoVector.push_back ( theHisto_g ->getTH1F() );
640  histoVector.push_back ( theHisto_dus->getTH1F() );
641  }
642  histoVector.push_back ( theHisto_c->getTH1F() );
643  histoVector.push_back ( theHisto_b->getTH1F() );
644  histoVector.push_back ( theHisto_ni->getTH1F() );
645  histoVector.push_back ( theHisto_dusg->getTH1F());
646  histoVector.push_back ( theHisto_pu->getTH1F());
647  }
648  return histoVector;
649 }
MonitorElement * theHisto_g
MonitorElement * theHisto_c
TH1F * getTH1F() const
MonitorElement * theHisto_u
MonitorElement * theHisto_all
MonitorElement * theHisto_pu
MonitorElement * theHisto_d
MonitorElement * theHisto_dusg
MonitorElement * theHisto_b
MonitorElement * theHisto_s
MonitorElement * theHisto_dus
MonitorElement * theHisto_ni
template<class T>
TH1F* FlavourHistograms< T >::histo_all ( ) const
inline

Definition at line 82 of file FlavourHistorgrams.h.

References MonitorElement::getTH1F(), and FlavourHistograms< T >::theHisto_all.

Referenced by FlavourHistograms< T >::divide().

82 { return theHisto_all->getTH1F() ; }
TH1F * getTH1F() const
MonitorElement * theHisto_all
template<class T>
TH1F* FlavourHistograms< T >::histo_b ( ) const
inline

Definition at line 87 of file FlavourHistorgrams.h.

References MonitorElement::getTH1F(), and FlavourHistograms< T >::theHisto_b.

Referenced by FlavourHistograms< T >::divide().

87 { return theHisto_b->getTH1F() ; }
TH1F * getTH1F() const
MonitorElement * theHisto_b
template<class T>
TH1F* FlavourHistograms< T >::histo_c ( ) const
inline

Definition at line 86 of file FlavourHistorgrams.h.

References MonitorElement::getTH1F(), and FlavourHistograms< T >::theHisto_c.

Referenced by FlavourHistograms< T >::divide().

86 { return theHisto_c->getTH1F() ; }
MonitorElement * theHisto_c
TH1F * getTH1F() const
template<class T>
TH1F* FlavourHistograms< T >::histo_d ( ) const
inline

Definition at line 83 of file FlavourHistorgrams.h.

References MonitorElement::getTH1F(), and FlavourHistograms< T >::theHisto_d.

Referenced by FlavourHistograms< T >::divide().

83 { return theHisto_d ->getTH1F() ; }
TH1F * getTH1F() const
MonitorElement * theHisto_d
template<class T>
TH1F* FlavourHistograms< T >::histo_dus ( ) const
inline

Definition at line 90 of file FlavourHistorgrams.h.

References MonitorElement::getTH1F(), and FlavourHistograms< T >::theHisto_dus.

Referenced by FlavourHistograms< T >::divide().

90 { return theHisto_dus->getTH1F() ; }
TH1F * getTH1F() const
MonitorElement * theHisto_dus
template<class T>
TH1F* FlavourHistograms< T >::histo_dusg ( ) const
inline

Definition at line 91 of file FlavourHistorgrams.h.

References MonitorElement::getTH1F(), and FlavourHistograms< T >::theHisto_dusg.

Referenced by FlavourHistograms< T >::divide().

91 { return theHisto_dusg->getTH1F() ; }
TH1F * getTH1F() const
MonitorElement * theHisto_dusg
template<class T>
TH1F* FlavourHistograms< T >::histo_g ( ) const
inline

Definition at line 88 of file FlavourHistorgrams.h.

References MonitorElement::getTH1F(), and FlavourHistograms< T >::theHisto_g.

Referenced by FlavourHistograms< T >::divide().

88 { return theHisto_g->getTH1F() ; }
MonitorElement * theHisto_g
TH1F * getTH1F() const
template<class T>
TH1F* FlavourHistograms< T >::histo_ni ( ) const
inline

Definition at line 89 of file FlavourHistorgrams.h.

References MonitorElement::getTH1F(), and FlavourHistograms< T >::theHisto_ni.

Referenced by FlavourHistograms< T >::divide().

89 { return theHisto_ni->getTH1F() ; }
TH1F * getTH1F() const
MonitorElement * theHisto_ni
template<class T>
TH1F* FlavourHistograms< T >::histo_pu ( ) const
inline
template<class T>
TH1F* FlavourHistograms< T >::histo_s ( ) const
inline

Definition at line 85 of file FlavourHistorgrams.h.

References MonitorElement::getTH1F(), and FlavourHistograms< T >::theHisto_s.

Referenced by FlavourHistograms< T >::divide().

85 { return theHisto_s->getTH1F() ; }
TH1F * getTH1F() const
MonitorElement * theHisto_s
template<class T>
TH1F* FlavourHistograms< T >::histo_u ( ) const
inline

Definition at line 84 of file FlavourHistorgrams.h.

References MonitorElement::getTH1F(), and FlavourHistograms< T >::theHisto_u.

Referenced by FlavourHistograms< T >::divide().

84 { return theHisto_u->getTH1F() ; }
TH1F * getTH1F() const
MonitorElement * theHisto_u
template<class T>
int FlavourHistograms< T >::indexToPlot ( ) const
inline
template<class T>
double FlavourHistograms< T >::lowerBound ( ) const
inline

Definition at line 71 of file FlavourHistorgrams.h.

References FlavourHistograms< T >::theLowerBound.

Referenced by EffPurFromHistos::EffPurFromHistos().

71 { return theLowerBound ; }
template<class T>
int FlavourHistograms< T >::maxDimension ( ) const
inline
template<class T>
int FlavourHistograms< T >::nBins ( ) const
inline

Definition at line 70 of file FlavourHistorgrams.h.

References FlavourHistograms< T >::theNBins.

Referenced by EffPurFromHistos::EffPurFromHistos().

70 { return theNBins ; }
template<class T >
void FlavourHistograms< T >::plot ( TPad *  theCanvas = 0)

Definition at line 367 of file FlavourHistorgrams.h.

References cuy::col, MonitorElement::getTH1F(), trackerHits::histo, mps_fire::i, RecoTauValidation_cfi::lineStyle, RecoTauValidation_cfi::lineWidth, RecoTauValidation_cfi::markerSize, RecoTauValidation_cfi::markerStyle, SiStripPI::max, nullptr, RecoBTag::setTDRStyle(), FlavourHistograms< T >::theBaseNameDescription, FlavourHistograms< T >::theHisto_b, FlavourHistograms< T >::theHisto_c, FlavourHistograms< T >::theHisto_dusg, FlavourHistograms< T >::theHisto_ni, FlavourHistograms< T >::theMax, FlavourHistograms< T >::theMin, FlavourHistograms< T >::thePlotFirst, FlavourHistograms< T >::thePlotLog, and FlavourHistograms< T >::thePlotNormalized.

Referenced by FlavourHistograms< T >::epsPlot().

367  {
368 
369 //fixme:
370  bool btppNI = false;
371  bool btppColour = true;
372 
373  if (theCanvas)
374  theCanvas->cd();
375 
376  RecoBTag::setTDRStyle()->cd();
377  gPad->UseCurrentStyle();
378 // if ( !btppTitle ) gStyle->SetOptTitle ( 0 ) ;
379 //
380 // // here: plot histograms in a canvas
381 // theCanvas = new TCanvas ( "C" + theBaseNameTitle , "C" + theBaseNameDescription ,
382 // btppXCanvas , btppYCanvas ) ;
383 // theCanvas->SetFillColor ( 0 ) ;
384 // theCanvas->cd ( 1 ) ;
385  gPad->SetLogy ( 0 ) ;
386  if ( thePlotLog ) gPad->SetLogy ( 1 ) ;
387  gPad->SetGridx ( 0 ) ;
388  gPad->SetGridy ( 0 ) ;
389  gPad->SetTitle ( nullptr ) ;
390 
391  MonitorElement * histo[4];
392  int col[4], lineStyle[4], markerStyle[4];
393  int lineWidth = 1 ;
394 
395  const double markerSize = gPad->GetWh() * gPad->GetHNDC() / 500.;
396 
397  // default (l)
398  histo[0] = theHisto_dusg ;
399  //CW histo_1 = theHisto_dus ;
400  histo[1] = theHisto_b ;
401  histo[2] = theHisto_c ;
402  histo[3]= nullptr ;
403 
404  double max = theMax;
405  if (theMax<=0.) {
406  max = theHisto_dusg->getTH1F()->GetMaximum();
407  if (theHisto_b->getTH1F()->GetMaximum() > max) max = theHisto_b->getTH1F()->GetMaximum();
408  if (theHisto_c->getTH1F()->GetMaximum() > max) max = theHisto_c->getTH1F()->GetMaximum();
409  }
410 
411  if (btppNI) {
412  histo[3] = theHisto_ni ;
413  if (theHisto_ni->getTH1F()->GetMaximum() > max) max = theHisto_ni->getTH1F()->GetMaximum();
414  }
415 
416  if ( btppColour ) { // print colours
417  col[0] = 4 ;
418  col[1] = 2 ;
419  col[2] = 6 ;
420  col[3] = 3 ;
421  lineStyle[0] = 1 ;
422  lineStyle[1] = 1 ;
423  lineStyle[2] = 1 ;
424  lineStyle[3] = 1 ;
425  markerStyle[0] = 20 ;
426  markerStyle[1] = 21 ;
427  markerStyle[2] = 22 ;
428  markerStyle[3] = 23 ;
429  lineWidth = 1 ;
430  }
431  else { // different marker/line styles
432  col[1] = 1 ;
433  col[2] = 1 ;
434  col[3] = 1 ;
435  col[0] = 1 ;
436  lineStyle[0] = 2 ;
437  lineStyle[1] = 1 ;
438  lineStyle[2] = 3 ;
439  lineStyle[3] = 4 ;
440  markerStyle[0] = 20 ;
441  markerStyle[1] = 21 ;
442  markerStyle[2] = 22 ;
443  markerStyle[3] = 23 ;
444  }
445 
446  // if changing order (NI stays always last)
447 
448  // c to plot first
449  if ( thePlotFirst == "c" ) {
450  histo[0] = theHisto_c ;
451  if ( btppColour ) col[0] = 6 ;
452  if ( !btppColour ) lineStyle[0] = 3 ;
453  histo[2] = theHisto_dusg ;
454  if ( btppColour ) col[2] = 4 ;
455  if ( !btppColour ) lineStyle[2] = 2 ;
456  }
457 
458  // b to plot first
459  if ( thePlotFirst == "b" ) {
460  histo[0] = theHisto_b ;
461  if ( btppColour ) col[0] = 2 ;
462  if ( !btppColour ) lineStyle[0] = 1 ;
463  histo[1] = theHisto_dusg ;
464  if ( btppColour ) col[1] = 4 ;
465  if ( !btppColour ) lineStyle[1] = 2 ;
466  }
467 
468 
469  histo[0] ->getTH1F()->GetXaxis()->SetTitle ( theBaseNameDescription.c_str() ) ;
470  histo[0] ->getTH1F()->GetYaxis()->SetTitle ( "Arbitrary Units" ) ;
471  histo[0] ->getTH1F()->GetYaxis()->SetTitleOffset(1.25) ;
472 
473  for (int i=0; i != 4; ++i) {
474  if (histo[i]== nullptr ) continue;
475  histo[i] ->getTH1F()->SetStats ( false ) ;
476  histo[i] ->getTH1F()->SetLineStyle ( lineStyle[i] ) ;
477  histo[i] ->getTH1F()->SetLineWidth ( lineWidth ) ;
478  histo[i] ->getTH1F()->SetLineColor ( col[i] ) ;
479  histo[i] ->getTH1F()->SetMarkerStyle ( markerStyle[i] ) ;
480  histo[i] ->getTH1F()->SetMarkerColor ( col[i] ) ;
481  histo[i] ->getTH1F()->SetMarkerSize ( markerSize ) ;
482  }
483 
484  if ( thePlotNormalized ) {
485  if (histo[0]->getTH1F()->GetEntries() != 0) {histo[0]->getTH1F() ->DrawNormalized() ;} else { histo[0]->getTH1F() ->SetMaximum(1.0); histo[0] ->getTH1F()->Draw() ;}
486  if (histo[1]->getTH1F()->GetEntries() != 0) histo[1] ->getTH1F()->DrawNormalized("Same") ;
487  if (histo[2]->getTH1F()->GetEntries() != 0) histo[2]->getTH1F() ->DrawNormalized("Same") ;
488  if ((histo[3] != nullptr) && (histo[3]->getTH1F()->GetEntries() != 0)) histo[3] ->getTH1F()->DrawNormalized("Same") ;
489  }
490  else {
491  histo[0]->getTH1F()->SetMaximum(max*1.05);
492  if (theMin!=-1.) histo[0]->getTH1F()->SetMinimum(theMin);
493  histo[0]->getTH1F()->Draw() ;
494  histo[1]->getTH1F()->Draw("Same") ;
495  histo[2]->getTH1F()->Draw("Same") ;
496  if ( histo[3] != nullptr ) histo[3]->getTH1F()->Draw("Same") ;
497  }
498 
499 }
MonitorElement * theHisto_c
TH1F * getTH1F() const
TStyle * setTDRStyle()
Definition: Tools.cc:346
#define nullptr
MonitorElement * theHisto_dusg
std::string theBaseNameDescription
MonitorElement * theHisto_b
col
Definition: cuy.py:1008
MonitorElement * theHisto_ni
template<class T>
std::string FlavourHistograms< T >::plotFirst ( ) const
inline

Definition at line 76 of file FlavourHistorgrams.h.

References FlavourHistograms< T >::thePlotFirst.

76 { return thePlotFirst ; }
template<class T>
bool FlavourHistograms< T >::plotLog ( ) const
inline

Definition at line 74 of file FlavourHistorgrams.h.

References FlavourHistograms< T >::thePlotLog.

74 { return thePlotLog ; }
template<class T>
bool FlavourHistograms< T >::plotNormalized ( ) const
inline

Definition at line 75 of file FlavourHistorgrams.h.

References FlavourHistograms< T >::thePlotNormalized.

75 { return thePlotNormalized ; }
template<class T >
void FlavourHistograms< T >::setEfficiencyFlag ( )

Definition at line 554 of file FlavourHistorgrams.h.

References FlavourHistograms< T >::mcPlots_, MonitorElement::setEfficiencyFlag(), FlavourHistograms< T >::theHisto_all, FlavourHistograms< T >::theHisto_b, FlavourHistograms< T >::theHisto_c, FlavourHistograms< T >::theHisto_d, FlavourHistograms< T >::theHisto_dus, FlavourHistograms< T >::theHisto_dusg, FlavourHistograms< T >::theHisto_g, FlavourHistograms< T >::theHisto_ni, FlavourHistograms< T >::theHisto_pu, FlavourHistograms< T >::theHisto_s, and FlavourHistograms< T >::theHisto_u.

554  {
556  if (mcPlots_) {
557  if (mcPlots_>2 ) {
563  }
569  }
570 
571 }
MonitorElement * theHisto_g
MonitorElement * theHisto_c
MonitorElement * theHisto_u
MonitorElement * theHisto_all
void setEfficiencyFlag()
MonitorElement * theHisto_pu
MonitorElement * theHisto_d
MonitorElement * theHisto_dusg
MonitorElement * theHisto_b
MonitorElement * theHisto_s
MonitorElement * theHisto_dus
MonitorElement * theHisto_ni
template<class T>
void FlavourHistograms< T >::SetMaximum ( const double &  max)
inline
template<class T>
void FlavourHistograms< T >::SetMinimum ( const double &  min)
inline

Definition at line 64 of file FlavourHistorgrams.h.

References min(), and FlavourHistograms< T >::theMin.

64 { theMin = min;}
T min(T a, T b)
Definition: MathUtil.h:58
template<class T >
void FlavourHistograms< T >::settitle ( const char *  title)

Definition at line 347 of file FlavourHistorgrams.h.

References FlavourHistograms< T >::mcPlots_, MonitorElement::setAxisTitle(), FlavourHistograms< T >::theHisto_all, FlavourHistograms< T >::theHisto_b, FlavourHistograms< T >::theHisto_c, FlavourHistograms< T >::theHisto_d, FlavourHistograms< T >::theHisto_dus, FlavourHistograms< T >::theHisto_dusg, FlavourHistograms< T >::theHisto_g, FlavourHistograms< T >::theHisto_ni, FlavourHistograms< T >::theHisto_pu, FlavourHistograms< T >::theHisto_s, and FlavourHistograms< T >::theHisto_u.

Referenced by TrackIPHistograms< T >::settitle().

347  {
349  if (mcPlots_) {
350  if (mcPlots_>2 ) {
356  }
362  }
363 }
MonitorElement * theHisto_g
MonitorElement * theHisto_c
MonitorElement * theHisto_u
MonitorElement * theHisto_all
MonitorElement * theHisto_pu
MonitorElement * theHisto_d
MonitorElement * theHisto_dusg
MonitorElement * theHisto_b
MonitorElement * theHisto_s
void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
MonitorElement * theHisto_dus
MonitorElement * theHisto_ni
template<class T>
bool FlavourHistograms< T >::statistics ( ) const
inline

Definition at line 73 of file FlavourHistorgrams.h.

References FlavourHistograms< T >::theStatistics.

73 { return theStatistics ; }
template<class T>
double FlavourHistograms< T >::upperBound ( ) const
inline

Definition at line 72 of file FlavourHistorgrams.h.

References FlavourHistograms< T >::theUpperBound.

Referenced by EffPurFromHistos::EffPurFromHistos().

72 { return theUpperBound ; }

Member Data Documentation

template<class T>
unsigned int FlavourHistograms< T >::mcPlots_
private
template<class T>
int* FlavourHistograms< T >::theArrayDimension
protected
template<class T>
std::string FlavourHistograms< T >::theBaseNameDescription
protected
template<class T>
std::string FlavourHistograms< T >::theBaseNameTitle
protected
template<class T>
TCanvas* FlavourHistograms< T >::theCanvas
protected
template<class T>
MonitorElement* FlavourHistograms< T >::theHisto_all
protected
template<class T>
MonitorElement* FlavourHistograms< T >::theHisto_b
protected
template<class T>
MonitorElement* FlavourHistograms< T >::theHisto_c
protected
template<class T>
MonitorElement* FlavourHistograms< T >::theHisto_d
protected
template<class T>
MonitorElement* FlavourHistograms< T >::theHisto_dus
protected
template<class T>
MonitorElement* FlavourHistograms< T >::theHisto_dusg
protected
template<class T>
MonitorElement* FlavourHistograms< T >::theHisto_g
protected
template<class T>
MonitorElement* FlavourHistograms< T >::theHisto_ni
protected
template<class T>
MonitorElement* FlavourHistograms< T >::theHisto_pu
protected
template<class T>
MonitorElement* FlavourHistograms< T >::theHisto_s
protected
template<class T>
MonitorElement* FlavourHistograms< T >::theHisto_u
protected
template<class T>
int FlavourHistograms< T >::theIndexToPlot
protected
template<class T>
double FlavourHistograms< T >::theLowerBound
protected
template<class T>
double FlavourHistograms< T >::theMax
protected
template<class T>
int FlavourHistograms< T >::theMaxDimension
protected
template<class T>
double FlavourHistograms< T >::theMin
protected
template<class T>
int FlavourHistograms< T >::theNBins
protected
template<class T>
std::string FlavourHistograms< T >::thePlotFirst
protected
template<class T>
bool FlavourHistograms< T >::thePlotLog
protected
template<class T>
bool FlavourHistograms< T >::thePlotNormalized
protected
template<class T>
bool FlavourHistograms< T >::theStatistics
protected
template<class T>
double FlavourHistograms< T >::theUpperBound
protected