CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/src/DQMOffline/PFTau/interface/PFJetMonitor.h

Go to the documentation of this file.
00001 #ifndef DQMOffline_PFTau_PFJetMonitor_h
00002 #define DQMOffline_PFTau_PFJetMonitor_h
00003 
00004 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00005 #include "DQMOffline/PFTau/interface/Benchmark.h"
00006 #include "DQMOffline/PFTau/interface/CandidateBenchmark.h"
00007 #include "DQMOffline/PFTau/interface/MatchCandidateBenchmark.h"
00008 
00009 #include "DataFormats/JetReco/interface/BasicJetCollection.h"
00010 
00011 #include <vector>
00012 class PFJetMonitor : public Benchmark {
00013 
00014  public:
00015 
00016   PFJetMonitor( float dRMax = 0.3,
00017                       bool matchCharge = true, 
00018                       Benchmark::Mode mode=Benchmark::DEFAULT) 
00019     : 
00020     Benchmark(mode), 
00021     candBench_(mode), matchCandBench_(mode), 
00022     dRMax_(dRMax), matchCharge_(matchCharge) {}
00023   
00024   virtual ~PFJetMonitor();
00025   
00027   void setParameters( const edm::ParameterSet& parameterSet);
00028   
00030   void setDirectory(TDirectory* dir);
00031 
00033   void setup();
00034   
00036   void setup(const edm::ParameterSet & parameterSet);
00037 
00039   template< class T, class C>
00040   void fill(const T& jetCollection,
00041             const C& matchedJetCollection, float& minVal, float& maxVal);
00042 
00043 
00044   void fillOne(const reco::Jet& jet,
00045                const reco::Jet& matchedJet);
00046 
00047  protected:
00048   CandidateBenchmark      candBench_;
00049   MatchCandidateBenchmark matchCandBench_;
00050 
00051   TH2F*  delta_frac_VS_frac_muon_;
00052   TH2F*  delta_frac_VS_frac_photon_;
00053   TH2F*  delta_frac_VS_frac_electron_;
00054   TH2F*  delta_frac_VS_frac_charged_hadron_;
00055   TH2F*  delta_frac_VS_frac_neutral_hadron_;
00056 
00057   float dRMax_;
00058   bool  matchCharge_;
00059   bool  createPFractionHistos_;
00060 
00061 };
00062 
00063 #include "DQMOffline/PFTau/interface/Matchers.h"
00064 template< class T, class C>
00065 void PFJetMonitor::fill(const T& jetCollection,
00066                         const C& matchedJetCollection, float& minVal, float& maxVal) {
00067   
00068 
00069   std::vector<int> matchIndices;
00070   PFB::match( jetCollection, matchedJetCollection, matchIndices, 
00071               matchCharge_, dRMax_ );
00072 
00073   for (unsigned int i = 0; i < (jetCollection).size(); i++) {
00074     const reco::Jet& jet = jetCollection[i];
00075 
00076     if( !isInRange(jet.pt(), jet.eta(), jet.phi() ) ) continue;
00077     
00078     int iMatch = matchIndices[i];
00079     assert(iMatch< static_cast<int>(matchedJetCollection.size()));
00080 
00081     if( iMatch!=-1 ) {
00082       const reco::Candidate& matchedJet = matchedJetCollection[ iMatch ];
00083       if(!isInRange(matchedJet.pt(),matchedJet.eta(),matchedJet.phi() ) ) continue;
00084       float ptRes = (jet.pt() - matchedJet.pt())/matchedJet.pt();
00085       
00086       if (ptRes > maxVal) maxVal = ptRes;
00087       if (ptRes < minVal) minVal = ptRes;
00088  
00089       candBench_.fillOne(jet);
00090       matchCandBench_.fillOne(jet, matchedJetCollection[ iMatch ]);
00091       if (createPFractionHistos_) fillOne(jet, matchedJetCollection[ iMatch ]);
00092     }
00093   }
00094 }
00095 #endif