CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_1/src/DQMOffline/RecoB/interface/TrackIPHistograms.h

Go to the documentation of this file.
00001 #ifndef DQMOffline_RecoB_TrackIPHistograms_h
00002 #define DQMOffline_RecoB_TrackIPHistograms_h
00003 
00004 #include <string>
00005 
00006 #include "DataFormats/TrackReco/interface/Track.h"
00007 #include "DQMOffline/RecoB/interface/FlavourHistorgrams.h"
00008 
00009 template<class T>
00010 class TrackIPHistograms : public FlavourHistograms<T>
00011 {
00012   public:
00013 
00014   TrackIPHistograms(const std::string& baseNameTitle_ , const std::string& baseNameDescription_,
00015                     const int& nBins_, const double& lowerBound_, const double& upperBound_,
00016                     const bool& statistics, const bool& plotLog_, const bool& plotNormalized_,
00017                     const std::string& plotFirst_, const bool& update, const std::string& folder, const unsigned int& mc, const bool& quality);
00018 
00019   virtual ~TrackIPHistograms(){};
00020 
00021   void fill(const int& flavour, const reco::TrackBase::TrackQuality& quality, const T& variable, const bool& hasTrack) const;
00022   void fill(const int& flavour, const reco::TrackBase::TrackQuality& quality, const T& variable, const bool& hasTrack, const T & w) const;
00023  
00024   void fill(const int& flavour, const reco::TrackBase::TrackQuality& quality, const T* variable, const bool& hasTrack) const;
00025   void fill(const int& flavour, const reco::TrackBase::TrackQuality& quality, const T* variable, const bool& hasTrack, const T & w) const;
00026 
00027   void settitle(const char* title);
00028 
00029   protected:
00030 
00031   void fillVariable ( const reco::TrackBase::TrackQuality& qual, const T & var, const bool& hasTrack) const;
00032   void fillVariable ( const reco::TrackBase::TrackQuality& qual, const T & var, const bool& hasTrack, const T & w) const;
00033 
00034   bool quality_;
00035 
00036   MonitorElement *theQual_undefined;
00037   MonitorElement *theQual_loose;
00038   MonitorElement *theQual_tight;
00039   MonitorElement *theQual_highpur;
00040 
00041   private:
00042   TrackIPHistograms(){}
00043 };
00044 
00045 template <class T>
00046 TrackIPHistograms<T>::TrackIPHistograms (const std::string& baseNameTitle_, const std::string& baseNameDescription_,
00047                                          const int& nBins_, const double& lowerBound_, const double& upperBound_,
00048                                          const bool& statistics_, const bool& plotLog_, const bool& plotNormalized_,
00049                                          const std::string& plotFirst_, const bool& update, const std::string& folder, const unsigned int& mc, const bool& quality) :
00050   FlavourHistograms<T>(baseNameTitle_, baseNameDescription_, nBins_, lowerBound_, upperBound_, statistics_, plotLog_, plotNormalized_,
00051                        plotFirst_, update, folder, mc), quality_(quality)
00052 {
00053   if(quality_) {
00054     if(!update) {
00055       HistoProviderDQM prov("Btag",folder);
00056       theQual_undefined = prov.book1D( baseNameTitle_ + "QualUnDef" , baseNameDescription_ + " Undefined Quality", nBins_, lowerBound_, upperBound_);
00057       theQual_loose = prov.book1D( baseNameTitle_ + "QualLoose" , baseNameDescription_ + " Loose Quality", nBins_, lowerBound_, upperBound_);
00058       theQual_tight = prov.book1D( baseNameTitle_ + "QualTight" , baseNameDescription_ + " Tight Quality", nBins_, lowerBound_, upperBound_);
00059       theQual_highpur = prov.book1D( baseNameTitle_ + "QualHighPur" , baseNameDescription_ + " High Purity Quality", nBins_, lowerBound_, upperBound_);
00060 
00061       if( statistics_ ) {
00062         theQual_undefined->getTH1F()->Sumw2();
00063         theQual_loose->getTH1F()->Sumw2();
00064         theQual_tight->getTH1F()->Sumw2();
00065         theQual_highpur->getTH1F()->Sumw2();
00066       }
00067     } else {
00068       HistoProviderDQM prov("Btag",folder);
00069       theQual_undefined = prov.access(baseNameTitle_ + "QualUnDef");
00070       theQual_loose = prov.access(baseNameTitle_ + "QualLoose");
00071       theQual_tight = prov.access(baseNameTitle_ + "QualTight");
00072       theQual_highpur = prov.access(baseNameTitle_ + "QualHighPur");
00073     }
00074   }
00075 }
00076 
00077 template <class T>
00078 void TrackIPHistograms<T>::fill(const int& flavour, const reco::TrackBase::TrackQuality& quality, const T& variable, const bool& hasTrack) const
00079 {
00080   FlavourHistograms<T>::fill(flavour, variable);
00081   if(quality_)
00082     fillVariable(quality, variable, hasTrack);
00083 }
00084 
00085 template <class T>
00086 void TrackIPHistograms<T>::fill(const int& flavour, const reco::TrackBase::TrackQuality& quality, const T& variable, const bool& hasTrack, const T & w) const
00087 {
00088   FlavourHistograms<T>::fill(flavour, variable , w);
00089   if(quality_)
00090     fillVariable(quality, variable, hasTrack, w);
00091 }
00092 
00093 template <class T>
00094 void TrackIPHistograms<T>::fill(const int& flavour, const reco::TrackBase::TrackQuality& quality, const T* variable, const bool& hasTrack) const
00095 {
00096   const int* theArrayDimension = FlavourHistograms<T>::arrayDimension();
00097   const int& theMaxDimension = FlavourHistograms<T>::maxDimension();
00098   const int& theIndexToPlot = FlavourHistograms<T>::indexToPlot();
00099 
00100   FlavourHistograms<T>::fill(flavour, variable);
00101   if( theArrayDimension == 0 && quality_) {
00102     fillVariable( quality, *variable);
00103   } else {
00104       int iMax = (*theArrayDimension > theMaxDimension) ? theMaxDimension : *theArrayDimension ;
00105       for(int i = 0; i != iMax; ++i) {
00106         if( quality_ && (( theIndexToPlot < 0) || ( i == theIndexToPlot)) ) {
00107           fillVariable ( flavour , *(variable + i), hasTrack);
00108         }
00109       }
00110 
00111       if(theIndexToPlot >= iMax && quality_) {
00112         const T& theZero = static_cast<T> (0.0);
00113         fillVariable ( quality, theZero, hasTrack);
00114       }
00115   }
00116 }
00117 
00118 template <class T>
00119 void TrackIPHistograms<T>::fill(const int& flavour, const reco::TrackBase::TrackQuality& quality, const T* variable, const bool& hasTrack, const T & w) const
00120 {
00121   const int* theArrayDimension = FlavourHistograms<T>::arrayDimension();
00122   const int& theMaxDimension = FlavourHistograms<T>::maxDimension();
00123   const int& theIndexToPlot = FlavourHistograms<T>::indexToPlot();
00124 
00125   FlavourHistograms<T>::fill(flavour, variable ,w);
00126   if( theArrayDimension == 0 && quality_) {
00127     fillVariable( quality, *variable,w);
00128   } else {
00129       int iMax = (*theArrayDimension > theMaxDimension) ? theMaxDimension : *theArrayDimension ;
00130       for(int i = 0; i != iMax; ++i) {
00131         if( quality_ && (( theIndexToPlot < 0) || ( i == theIndexToPlot)) ) {
00132           fillVariable ( flavour , *(variable + i), hasTrack,w);
00133         }
00134       }
00135 
00136       if(theIndexToPlot >= iMax && quality_) {
00137         const T& theZero = static_cast<T> (0.0);
00138         fillVariable ( quality, theZero, hasTrack,w);
00139       }
00140   }
00141 }
00142 
00143 template <class T>
00144 void TrackIPHistograms<T>::settitle(const char* title)
00145 {
00146   FlavourHistograms<T>::settitle(title);
00147   theQual_undefined->setAxisTitle(title);
00148   theQual_loose->setAxisTitle(title);
00149   theQual_tight->setAxisTitle(title);
00150   theQual_highpur->setAxisTitle(title);
00151 }
00152 
00153 template<class T>
00154 void TrackIPHistograms<T>::fillVariable( const reco::TrackBase::TrackQuality& qual, const T& var, const bool& hasTrack) const
00155 {
00156   if(!hasTrack || !quality_) return;
00157 
00158   switch(qual) {
00159     case reco::TrackBase::loose:
00160       theQual_loose->Fill(var);
00161       break;
00162     case reco::TrackBase::tight:
00163       theQual_tight->Fill(var);
00164       break;
00165     case reco::TrackBase::highPurity:
00166       theQual_highpur->Fill(var);
00167       break;
00168     default:
00169       theQual_undefined->Fill(var);
00170       break;
00171   }
00172 }
00173 
00174 template<class T>
00175 void TrackIPHistograms<T>::fillVariable( const reco::TrackBase::TrackQuality& qual, const T& var, const bool& hasTrack, const T & w) const
00176 {
00177   if(!hasTrack || !quality_) return;
00178 
00179   switch(qual) {
00180     case reco::TrackBase::loose:
00181       theQual_loose->Fill(var,w);
00182       break;
00183     case reco::TrackBase::tight:
00184       theQual_tight->Fill(var,w);
00185       break;
00186     case reco::TrackBase::highPurity:
00187       theQual_highpur->Fill(var,w);
00188       break;
00189     default:
00190       theQual_undefined->Fill(var,w);
00191       break;
00192   }
00193 }
00194 
00195 #endif