CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions
TrackIPHistograms< T > Class Template Reference

#include <TrackIPHistograms.h>

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

Public Member Functions

void fill (const int &flavour, const reco::TrackBase::TrackQuality &quality, const T &variable, const bool &hasTrack) const
 
void fill (const int &flavour, const reco::TrackBase::TrackQuality &quality, const T &variable, const bool &hasTrack, const T &w) const
 
void fill (const int &flavour, const reco::TrackBase::TrackQuality &quality, const T *variable, const bool &hasTrack) const
 
void fill (const int &flavour, const reco::TrackBase::TrackQuality &quality, const T *variable, const bool &hasTrack, const T &w) const
 
void settitle (const char *title)
 
 TrackIPHistograms (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, const bool &quality, DQMStore::IGetter &iget)
 
 TrackIPHistograms (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, const bool &quality, DQMStore::IBooker &ibook)
 
virtual ~TrackIPHistograms ()
 
- Public Member Functions inherited from FlavourHistograms< T >
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

void fillVariable (const reco::TrackBase::TrackQuality &qual, const T &var, const bool &hasTrack) const
 
void fillVariable (const reco::TrackBase::TrackQuality &qual, const T &var, const bool &hasTrack, const T &w) const
 
- Protected Member Functions inherited from FlavourHistograms< T >
double ClopperPearsonUnc (TH1F *num, TH1F *den, int bin)
 
void ComputeEfficiency (TH1F *num, TH1F *den, int bin)
 
void fillVariable (const int &flavour, const T &var, const T &w) const
 

Protected Attributes

bool quality_
 
MonitorElementtheQual_highpur
 
MonitorElementtheQual_loose
 
MonitorElementtheQual_tight
 
MonitorElementtheQual_undefined
 
- Protected Attributes inherited from FlavourHistograms< T >
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

 TrackIPHistograms ()
 

Detailed Description

template<class T>
class TrackIPHistograms< T >

Definition at line 10 of file TrackIPHistograms.h.

Constructor & Destructor Documentation

template<class T >
TrackIPHistograms< T >::TrackIPHistograms ( 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,
const bool &  quality,
DQMStore::IGetter iget 
)

Definition at line 76 of file TrackIPHistograms.h.

References DQMStore::IGetter::get(), TrackIPHistograms< T >::quality_, TrackIPHistograms< T >::theQual_highpur, TrackIPHistograms< T >::theQual_loose, TrackIPHistograms< T >::theQual_tight, and TrackIPHistograms< T >::theQual_undefined.

79  :
80 FlavourHistograms<T>(baseNameTitle_, baseNameDescription_, nBins_, lowerBound_, upperBound_,
81  plotFirst_, folder, mc, iget), quality_(quality)
82 {
83  if(quality_) {
84  theQual_undefined = iget.get("Btag/" + folder + "/" + baseNameTitle_ + "QualUnDef");
85  theQual_loose = iget.get("Btag/" + folder + "/" + baseNameTitle_ + "QualLoose");
86  theQual_tight = iget.get("Btag/" + folder + "/" + baseNameTitle_ + "QualTight");
87  theQual_highpur = iget.get("Btag/" + folder + "/" + baseNameTitle_ + "QualHighPur");
88  }
89 }
MonitorElement * theQual_tight
MonitorElement * get(const std::string &path)
Definition: DQMStore.cc:302
MonitorElement * theQual_loose
MonitorElement * theQual_highpur
MonitorElement * theQual_undefined
template<class T >
TrackIPHistograms< T >::TrackIPHistograms ( 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,
const bool &  quality,
DQMStore::IBooker ibook 
)

Definition at line 52 of file TrackIPHistograms.h.

References HistoProviderDQM::book1D(), MonitorElement::getTH1F(), TrackIPHistograms< T >::quality_, TrackIPHistograms< T >::theQual_highpur, TrackIPHistograms< T >::theQual_loose, TrackIPHistograms< T >::theQual_tight, and TrackIPHistograms< T >::theQual_undefined.

56  :
57  FlavourHistograms<T>(baseNameTitle_, baseNameDescription_, nBins_, lowerBound_, upperBound_, statistics_, plotLog_, plotNormalized_,
58  plotFirst_, folder, mc, ibook), quality_(quality)
59 {
60  if(quality_) {
61  HistoProviderDQM prov("Btag",folder,ibook);
62  theQual_undefined = prov.book1D( baseNameTitle_ + "QualUnDef" , baseNameDescription_ + " Undefined Quality", nBins_, lowerBound_, upperBound_);
63  theQual_loose = prov.book1D( baseNameTitle_ + "QualLoose" , baseNameDescription_ + " Loose Quality", nBins_, lowerBound_, upperBound_);
64  theQual_tight = prov.book1D( baseNameTitle_ + "QualTight" , baseNameDescription_ + " Tight Quality", nBins_, lowerBound_, upperBound_);
65  theQual_highpur = prov.book1D( baseNameTitle_ + "QualHighPur" , baseNameDescription_ + " High Purity Quality", nBins_, lowerBound_, upperBound_);
66 
67  if( statistics_ ) {
68  theQual_undefined->getTH1F()->Sumw2();
69  theQual_loose->getTH1F()->Sumw2();
70  theQual_tight->getTH1F()->Sumw2();
71  theQual_highpur->getTH1F()->Sumw2();
72  }
73  }
74 }
MonitorElement * theQual_tight
MonitorElement * theQual_loose
MonitorElement * theQual_highpur
TH1F * getTH1F(void) const
MonitorElement * theQual_undefined
template<class T>
virtual TrackIPHistograms< T >::~TrackIPHistograms ( )
inlinevirtual

Definition at line 25 of file TrackIPHistograms.h.

25 {};
template<class T>
TrackIPHistograms< T >::TrackIPHistograms ( )
inlineprivate

Definition at line 48 of file TrackIPHistograms.h.

48 {}

Member Function Documentation

template<class T>
void TrackIPHistograms< T >::fill ( const int &  flavour,
const reco::TrackBase::TrackQuality quality,
const T variable,
const bool &  hasTrack 
) const

Definition at line 92 of file TrackIPHistograms.h.

References FlavourHistograms< T >::fill().

93 {
95  if(quality_)
96  fillVariable(quality, variable, hasTrack);
97 }
void fill(const int &flavour, const T &variable) const
int flavour(const Candidate &part)
Definition: pdgIdUtils.h:31
void fillVariable(const reco::TrackBase::TrackQuality &qual, const T &var, const bool &hasTrack) const
template<class T>
void TrackIPHistograms< T >::fill ( const int &  flavour,
const reco::TrackBase::TrackQuality quality,
const T variable,
const bool &  hasTrack,
const T w 
) const

Definition at line 100 of file TrackIPHistograms.h.

References FlavourHistograms< T >::fill().

101 {
102  FlavourHistograms<T>::fill(flavour, variable , w);
103  if(quality_)
104  fillVariable(quality, variable, hasTrack, w);
105 }
const double w
Definition: UKUtility.cc:23
void fill(const int &flavour, const T &variable) const
int flavour(const Candidate &part)
Definition: pdgIdUtils.h:31
void fillVariable(const reco::TrackBase::TrackQuality &qual, const T &var, const bool &hasTrack) const
template<class T>
void TrackIPHistograms< T >::fill ( const int &  flavour,
const reco::TrackBase::TrackQuality quality,
const T variable,
const bool &  hasTrack 
) const

Definition at line 108 of file TrackIPHistograms.h.

References FlavourHistograms< T >::arrayDimension(), FlavourHistograms< T >::fill(), i, FlavourHistograms< T >::indexToPlot(), and FlavourHistograms< T >::maxDimension().

109 {
113 
115  if( theArrayDimension == 0 && quality_) {
116  fillVariable( quality, *variable);
117  } else {
118  int iMax = (*theArrayDimension > theMaxDimension) ? theMaxDimension : *theArrayDimension ;
119  for(int i = 0; i != iMax; ++i) {
120  if( quality_ && (( theIndexToPlot < 0) || ( i == theIndexToPlot)) ) {
121  fillVariable ( flavour , *(variable + i), hasTrack);
122  }
123  }
124 
125  if(theIndexToPlot >= iMax && quality_) {
126  const T& theZero = static_cast<T> (0.0);
127  fillVariable ( quality, theZero, hasTrack);
128  }
129  }
130 }
int i
Definition: DBlmapReader.cc:9
int maxDimension() const
void fill(const int &flavour, const T &variable) const
int * arrayDimension() const
int indexToPlot() const
long double T
int flavour(const Candidate &part)
Definition: pdgIdUtils.h:31
void fillVariable(const reco::TrackBase::TrackQuality &qual, const T &var, const bool &hasTrack) const
template<class T>
void TrackIPHistograms< T >::fill ( const int &  flavour,
const reco::TrackBase::TrackQuality quality,
const T variable,
const bool &  hasTrack,
const T w 
) const

Definition at line 133 of file TrackIPHistograms.h.

References FlavourHistograms< T >::arrayDimension(), FlavourHistograms< T >::fill(), i, FlavourHistograms< T >::indexToPlot(), and FlavourHistograms< T >::maxDimension().

134 {
138 
140  if( theArrayDimension == 0 && quality_) {
141  fillVariable( quality, *variable,w);
142  } else {
143  int iMax = (*theArrayDimension > theMaxDimension) ? theMaxDimension : *theArrayDimension ;
144  for(int i = 0; i != iMax; ++i) {
145  if( quality_ && (( theIndexToPlot < 0) || ( i == theIndexToPlot)) ) {
146  fillVariable ( flavour , *(variable + i), hasTrack,w);
147  }
148  }
149 
150  if(theIndexToPlot >= iMax && quality_) {
151  const T& theZero = static_cast<T> (0.0);
152  fillVariable ( quality, theZero, hasTrack,w);
153  }
154  }
155 }
int i
Definition: DBlmapReader.cc:9
const double w
Definition: UKUtility.cc:23
int maxDimension() const
void fill(const int &flavour, const T &variable) const
int * arrayDimension() const
int indexToPlot() const
long double T
int flavour(const Candidate &part)
Definition: pdgIdUtils.h:31
void fillVariable(const reco::TrackBase::TrackQuality &qual, const T &var, const bool &hasTrack) const
template<class T>
void TrackIPHistograms< T >::fillVariable ( const reco::TrackBase::TrackQuality qual,
const T var,
const bool &  hasTrack 
) const
protected

Definition at line 168 of file TrackIPHistograms.h.

References reco::TrackBase::highPurity, reco::TrackBase::loose, and reco::TrackBase::tight.

169 {
170  if(!hasTrack || !quality_) return;
171 
172  switch(qual) {
175  break;
178  break;
181  break;
182  default:
184  break;
185  }
186 }
MonitorElement * theQual_tight
void Fill(long long x)
MonitorElement * theQual_loose
MonitorElement * theQual_highpur
MonitorElement * theQual_undefined
template<class T>
void TrackIPHistograms< T >::fillVariable ( const reco::TrackBase::TrackQuality qual,
const T var,
const bool &  hasTrack,
const T w 
) const
protected

Definition at line 189 of file TrackIPHistograms.h.

References reco::TrackBase::highPurity, reco::TrackBase::loose, and reco::TrackBase::tight.

190 {
191  if(!hasTrack || !quality_) return;
192 
193  switch(qual) {
196  break;
199  break;
202  break;
203  default:
205  break;
206  }
207 }
MonitorElement * theQual_tight
const double w
Definition: UKUtility.cc:23
void Fill(long long x)
MonitorElement * theQual_loose
MonitorElement * theQual_highpur
MonitorElement * theQual_undefined
template<class T >
void TrackIPHistograms< T >::settitle ( const char *  title)

Definition at line 158 of file TrackIPHistograms.h.

References FlavourHistograms< T >::settitle().

159 {
165 }
MonitorElement * theQual_tight
MonitorElement * theQual_loose
void settitle(const char *title)
MonitorElement * theQual_highpur
MonitorElement * theQual_undefined
void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)

Member Data Documentation

template<class T>
bool TrackIPHistograms< T >::quality_
protected

Definition at line 40 of file TrackIPHistograms.h.

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

template<class T>
MonitorElement* TrackIPHistograms< T >::theQual_highpur
protected

Definition at line 45 of file TrackIPHistograms.h.

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

template<class T>
MonitorElement* TrackIPHistograms< T >::theQual_loose
protected

Definition at line 43 of file TrackIPHistograms.h.

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

template<class T>
MonitorElement* TrackIPHistograms< T >::theQual_tight
protected

Definition at line 44 of file TrackIPHistograms.h.

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

template<class T>
MonitorElement* TrackIPHistograms< T >::theQual_undefined
protected

Definition at line 42 of file TrackIPHistograms.h.

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