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 bool &statistics, const bool &plotLog_, const bool &plotNormalized_, const std::string &plotFirst_, const bool &update, const std::string &folder, const unsigned int &mc, const bool &quality)
 
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) const
 
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 bool &statistics_, const bool &plotLog_, const bool &plotNormalized_, const std::string &plotFirst_, const bool &update, const std::string &folder, const unsigned int &mc)
 
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_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 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 >
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_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 bool &  statistics,
const bool &  plotLog_,
const bool &  plotNormalized_,
const std::string &  plotFirst_,
const bool &  update,
const std::string &  folder,
const unsigned int &  mc,
const bool &  quality 
)

Definition at line 46 of file TrackIPHistograms.h.

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

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

Definition at line 19 of file TrackIPHistograms.h.

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

Definition at line 42 of file TrackIPHistograms.h.

42 {}

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 78 of file TrackIPHistograms.h.

References FlavourHistograms< T >::fill().

Referenced by TrackIPTagPlotter::analyzeTag().

79 {
81  if(quality_)
82  fillVariable(quality, variable, hasTrack);
83 }
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 86 of file TrackIPHistograms.h.

References FlavourHistograms< T >::fill().

87 {
89  if(quality_)
90  fillVariable(quality, variable, hasTrack, w);
91 }
void fill(const int &flavour, const T &variable) const
T w() 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 94 of file TrackIPHistograms.h.

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

95 {
99 
101  if( theArrayDimension == 0 && quality_) {
102  fillVariable( quality, *variable);
103  } else {
104  int iMax = (*theArrayDimension > theMaxDimension) ? theMaxDimension : *theArrayDimension ;
105  for(int i = 0; i != iMax; ++i) {
106  if( quality_ && (( theIndexToPlot < 0) || ( i == theIndexToPlot)) ) {
107  fillVariable ( flavour , *(variable + i), hasTrack);
108  }
109  }
110 
111  if(theIndexToPlot >= iMax && quality_) {
112  const T& theZero = static_cast<T> (0.0);
113  fillVariable ( quality, theZero, hasTrack);
114  }
115  }
116 }
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 119 of file TrackIPHistograms.h.

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

120 {
124 
126  if( theArrayDimension == 0 && quality_) {
127  fillVariable( quality, *variable,w);
128  } else {
129  int iMax = (*theArrayDimension > theMaxDimension) ? theMaxDimension : *theArrayDimension ;
130  for(int i = 0; i != iMax; ++i) {
131  if( quality_ && (( theIndexToPlot < 0) || ( i == theIndexToPlot)) ) {
132  fillVariable ( flavour , *(variable + i), hasTrack,w);
133  }
134  }
135 
136  if(theIndexToPlot >= iMax && quality_) {
137  const T& theZero = static_cast<T> (0.0);
138  fillVariable ( quality, theZero, hasTrack,w);
139  }
140  }
141 }
int i
Definition: DBlmapReader.cc:9
int maxDimension() const
void fill(const int &flavour, const T &variable) const
int * arrayDimension() const
int indexToPlot() const
T w() 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 154 of file TrackIPHistograms.h.

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

155 {
156  if(!hasTrack || !quality_) return;
157 
158  switch(qual) {
160  theQual_loose->Fill(var);
161  break;
163  theQual_tight->Fill(var);
164  break;
166  theQual_highpur->Fill(var);
167  break;
168  default:
169  theQual_undefined->Fill(var);
170  break;
171  }
172 }
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 175 of file TrackIPHistograms.h.

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

176 {
177  if(!hasTrack || !quality_) return;
178 
179  switch(qual) {
181  theQual_loose->Fill(var,w);
182  break;
184  theQual_tight->Fill(var,w);
185  break;
187  theQual_highpur->Fill(var,w);
188  break;
189  default:
190  theQual_undefined->Fill(var,w);
191  break;
192  }
193 }
MonitorElement * theQual_tight
void Fill(long long x)
MonitorElement * theQual_loose
MonitorElement * theQual_highpur
MonitorElement * theQual_undefined
T w() const
template<class T >
void TrackIPHistograms< T >::settitle ( const char *  title)

Definition at line 144 of file TrackIPHistograms.h.

References FlavourHistograms< T >::settitle().

145 {
151 }
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 34 of file TrackIPHistograms.h.

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

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

Definition at line 39 of file TrackIPHistograms.h.

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

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

Definition at line 37 of file TrackIPHistograms.h.

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

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

Definition at line 38 of file TrackIPHistograms.h.

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

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

Definition at line 36 of file TrackIPHistograms.h.

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