Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-09-11 04:32:51

0001 #ifndef DQMOffline_RecoB_TrackIPHistograms_h
0002 #define DQMOffline_RecoB_TrackIPHistograms_h
0003 
0004 #include <string>
0005 
0006 #include "DataFormats/TrackReco/interface/Track.h"
0007 #include "DQMOffline/RecoB/interface/FlavourHistorgrams.h"
0008 
0009 template <class T>
0010 class TrackIPHistograms : public FlavourHistograms<T> {
0011 public:
0012   typedef dqm::legacy::DQMStore DQMStore;
0013   typedef dqm::legacy::MonitorElement MonitorElement;
0014 
0015   TrackIPHistograms(const std::string& baseNameTitle_,
0016                     const std::string& baseNameDescription_,
0017                     const int& nBins_,
0018                     const double& lowerBound_,
0019                     const double& upperBound_,
0020                     const std::string& plotFirst_,
0021                     const std::string& folder,
0022                     const unsigned int& mc,
0023                     const bool& quality,
0024                     DQMStore::IGetter& iget);
0025 
0026   TrackIPHistograms(const std::string& baseNameTitle_,
0027                     const std::string& baseNameDescription_,
0028                     const int& nBins_,
0029                     const double& lowerBound_,
0030                     const double& upperBound_,
0031                     const bool& statistics,
0032                     const bool& plotLog_,
0033                     const bool& plotNormalized_,
0034                     const std::string& plotFirst_,
0035                     const std::string& folder,
0036                     const unsigned int& mc,
0037                     const bool& quality,
0038                     DQMStore::IBooker& ibook);
0039 
0040   ~TrackIPHistograms() override {}
0041 
0042   void fill(const int& flavour,
0043             const reco::TrackBase::TrackQuality& quality,
0044             const T& variable,
0045             const bool& hasTrack) const;
0046   void fill(const int& flavour,
0047             const reco::TrackBase::TrackQuality& quality,
0048             const T& variable,
0049             const bool& hasTrack,
0050             const T& w) const;
0051 
0052   void fill(const int& flavour,
0053             const reco::TrackBase::TrackQuality& quality,
0054             const T* variable,
0055             const bool& hasTrack) const;
0056   void fill(const int& flavour,
0057             const reco::TrackBase::TrackQuality& quality,
0058             const T* variable,
0059             const bool& hasTrack,
0060             const T& w) const;
0061 
0062   void settitle(const char* title);
0063 
0064 protected:
0065   void fillVariable(const reco::TrackBase::TrackQuality& qual, const T& var, const bool& hasTrack) const;
0066   void fillVariable(const reco::TrackBase::TrackQuality& qual, const T& var, const bool& hasTrack, const T& w) const;
0067 
0068   bool quality_;
0069 
0070   MonitorElement* theQual_undefined;
0071   MonitorElement* theQual_loose;
0072   MonitorElement* theQual_tight;
0073   MonitorElement* theQual_highpur;
0074 
0075 private:
0076   TrackIPHistograms() {}
0077 };
0078 
0079 template <class T>
0080 TrackIPHistograms<T>::TrackIPHistograms(const std::string& baseNameTitle_,
0081                                         const std::string& baseNameDescription_,
0082                                         const int& nBins_,
0083                                         const double& lowerBound_,
0084                                         const double& upperBound_,
0085                                         const bool& statistics_,
0086                                         const bool& plotLog_,
0087                                         const bool& plotNormalized_,
0088                                         const std::string& plotFirst_,
0089                                         const std::string& folder,
0090                                         const unsigned int& mc,
0091                                         const bool& quality,
0092                                         DQMStore::IBooker& ibook)
0093     : FlavourHistograms<T>(baseNameTitle_,
0094                            baseNameDescription_,
0095                            nBins_,
0096                            lowerBound_,
0097                            upperBound_,
0098                            statistics_,
0099                            plotLog_,
0100                            plotNormalized_,
0101                            plotFirst_,
0102                            folder,
0103                            mc,
0104                            ibook),
0105       quality_(quality) {
0106   if (quality_) {
0107     HistoProviderDQM prov("Btag", folder, ibook);
0108     theQual_undefined = prov.book1D(
0109         baseNameTitle_ + "QualUnDef", baseNameDescription_ + " Undefined Quality", nBins_, lowerBound_, upperBound_);
0110     theQual_loose = prov.book1D(
0111         baseNameTitle_ + "QualLoose", baseNameDescription_ + " Loose Quality", nBins_, lowerBound_, upperBound_);
0112     theQual_tight = prov.book1D(
0113         baseNameTitle_ + "QualTight", baseNameDescription_ + " Tight Quality", nBins_, lowerBound_, upperBound_);
0114     theQual_highpur = prov.book1D(
0115         baseNameTitle_ + "QualHighPur", baseNameDescription_ + " High Purity Quality", nBins_, lowerBound_, upperBound_);
0116 
0117     if (statistics_) {
0118       theQual_undefined->enableSumw2();
0119       theQual_loose->enableSumw2();
0120       theQual_tight->enableSumw2();
0121       theQual_highpur->enableSumw2();
0122     }
0123   }
0124 }
0125 template <class T>
0126 TrackIPHistograms<T>::TrackIPHistograms(const std::string& baseNameTitle_,
0127                                         const std::string& baseNameDescription_,
0128                                         const int& nBins_,
0129                                         const double& lowerBound_,
0130                                         const double& upperBound_,
0131                                         const std::string& plotFirst_,
0132                                         const std::string& folder,
0133                                         const unsigned int& mc,
0134                                         const bool& quality,
0135                                         DQMStore::IGetter& iget)
0136     : FlavourHistograms<T>(
0137           baseNameTitle_, baseNameDescription_, nBins_, lowerBound_, upperBound_, plotFirst_, folder, mc, iget),
0138       quality_(quality) {
0139   if (quality_) {
0140     theQual_undefined = iget.get("Btag/" + folder + "/" + baseNameTitle_ + "QualUnDef");
0141     theQual_loose = iget.get("Btag/" + folder + "/" + baseNameTitle_ + "QualLoose");
0142     theQual_tight = iget.get("Btag/" + folder + "/" + baseNameTitle_ + "QualTight");
0143     theQual_highpur = iget.get("Btag/" + folder + "/" + baseNameTitle_ + "QualHighPur");
0144   }
0145 }
0146 
0147 template <class T>
0148 void TrackIPHistograms<T>::fill(const int& flavour,
0149                                 const reco::TrackBase::TrackQuality& quality,
0150                                 const T& variable,
0151                                 const bool& hasTrack) const {
0152   FlavourHistograms<T>::fill(flavour, variable);
0153   if (quality_)
0154     fillVariable(quality, variable, hasTrack);
0155 }
0156 
0157 template <class T>
0158 void TrackIPHistograms<T>::fill(const int& flavour,
0159                                 const reco::TrackBase::TrackQuality& quality,
0160                                 const T& variable,
0161                                 const bool& hasTrack,
0162                                 const T& w) const {
0163   FlavourHistograms<T>::fill(flavour, variable, w);
0164   if (quality_)
0165     fillVariable(quality, variable, hasTrack, w);
0166 }
0167 
0168 template <class T>
0169 void TrackIPHistograms<T>::fill(const int& flavour,
0170                                 const reco::TrackBase::TrackQuality& quality,
0171                                 const T* variable,
0172                                 const bool& hasTrack) const {
0173   const int* theArrayDimension = FlavourHistograms<T>::arrayDimension();
0174   const int& theMaxDimension = FlavourHistograms<T>::maxDimension();
0175   const int& theIndexToPlot = FlavourHistograms<T>::indexToPlot();
0176 
0177   FlavourHistograms<T>::fill(flavour, variable);
0178   if (theArrayDimension == nullptr && quality_) {
0179     fillVariable(quality, *variable);
0180   } else {
0181     int iMax = (*theArrayDimension > theMaxDimension) ? theMaxDimension : *theArrayDimension;
0182     for (int i = 0; i != iMax; ++i) {
0183       if (quality_ && ((theIndexToPlot < 0) || (i == theIndexToPlot))) {
0184         fillVariable(flavour, *(variable + i), hasTrack);
0185       }
0186     }
0187 
0188     if (theIndexToPlot >= iMax && quality_) {
0189       const T& theZero = static_cast<T>(0.0);
0190       fillVariable(quality, theZero, hasTrack);
0191     }
0192   }
0193 }
0194 
0195 template <class T>
0196 void TrackIPHistograms<T>::fill(const int& flavour,
0197                                 const reco::TrackBase::TrackQuality& quality,
0198                                 const T* variable,
0199                                 const bool& hasTrack,
0200                                 const T& w) const {
0201   const int* theArrayDimension = FlavourHistograms<T>::arrayDimension();
0202   const int& theMaxDimension = FlavourHistograms<T>::maxDimension();
0203   const int& theIndexToPlot = FlavourHistograms<T>::indexToPlot();
0204 
0205   FlavourHistograms<T>::fill(flavour, variable, w);
0206   if (theArrayDimension == nullptr && quality_) {
0207     fillVariable(quality, *variable, w);
0208   } else {
0209     int iMax = (*theArrayDimension > theMaxDimension) ? theMaxDimension : *theArrayDimension;
0210     for (int i = 0; i != iMax; ++i) {
0211       if (quality_ && ((theIndexToPlot < 0) || (i == theIndexToPlot))) {
0212         fillVariable(flavour, *(variable + i), hasTrack, w);
0213       }
0214     }
0215 
0216     if (theIndexToPlot >= iMax && quality_) {
0217       const T& theZero = static_cast<T>(0.0);
0218       fillVariable(quality, theZero, hasTrack, w);
0219     }
0220   }
0221 }
0222 
0223 template <class T>
0224 void TrackIPHistograms<T>::settitle(const char* title) {
0225   FlavourHistograms<T>::settitle(title);
0226   theQual_undefined->setAxisTitle(title);
0227   theQual_loose->setAxisTitle(title);
0228   theQual_tight->setAxisTitle(title);
0229   theQual_highpur->setAxisTitle(title);
0230 }
0231 
0232 template <class T>
0233 void TrackIPHistograms<T>::fillVariable(const reco::TrackBase::TrackQuality& qual,
0234                                         const T& var,
0235                                         const bool& hasTrack) const {
0236   if (!hasTrack || !quality_)
0237     return;
0238 
0239   switch (qual) {
0240     case reco::TrackBase::loose:
0241       theQual_loose->Fill(var);
0242       break;
0243     case reco::TrackBase::tight:
0244       theQual_tight->Fill(var);
0245       break;
0246     case reco::TrackBase::highPurity:
0247       theQual_highpur->Fill(var);
0248       break;
0249     default:
0250       theQual_undefined->Fill(var);
0251       break;
0252   }
0253 }
0254 
0255 template <class T>
0256 void TrackIPHistograms<T>::fillVariable(const reco::TrackBase::TrackQuality& qual,
0257                                         const T& var,
0258                                         const bool& hasTrack,
0259                                         const T& w) const {
0260   if (!hasTrack || !quality_)
0261     return;
0262 
0263   switch (qual) {
0264     case reco::TrackBase::loose:
0265       theQual_loose->Fill(var, w);
0266       break;
0267     case reco::TrackBase::tight:
0268       theQual_tight->Fill(var, w);
0269       break;
0270     case reco::TrackBase::highPurity:
0271       theQual_highpur->Fill(var, w);
0272       break;
0273     default:
0274       theQual_undefined->Fill(var, w);
0275       break;
0276   }
0277 }
0278 
0279 #endif