CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch13/src/PhysicsTools/PatExamples/plugins/PatJetAnalyzer.cc

Go to the documentation of this file.
00001 #include <map>
00002 #include <string>
00003 #include <iomanip>
00004 #include <sstream>
00005 #include <iostream>
00006 
00007 #include "TH1.h"
00008 
00009 #include "FWCore/Framework/interface/Event.h"
00010 #include "FWCore/Framework/interface/EDAnalyzer.h"
00011 #include "FWCore/Utilities/interface/InputTag.h"
00012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00013 #include "FWCore/ServiceRegistry/interface/Service.h"
00014 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00015 
00016 #include "DataFormats/PatCandidates/interface/Jet.h"
00017 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00018 
00019 
00022 static const unsigned int MAXBIN=8;
00026 static const float BINS[]={30., 40., 50., 60., 70., 80., 100., 125., 150.};
00027 
00043 class PatJetAnalyzer : public edm::EDAnalyzer {
00044 
00045 public:
00047   explicit PatJetAnalyzer(const edm::ParameterSet& cfg);
00049   ~PatJetAnalyzer(){};
00050   
00051 private:
00053   virtual void analyze(const edm::Event& event, const edm::EventSetup& setup);
00054 
00056   bool booked(const std::string histName) const { return hists_.find(histName.c_str())!=hists_.end(); };
00058   void fill(const std::string histName, double value) const { if(booked(histName.c_str())) hists_.find(histName.c_str())->second->Fill(value); };
00059   // print jet pt for each level of energy corrections
00060   void print(edm::View<pat::Jet>::const_iterator& jet, unsigned int idx);
00061 
00062 private:  
00064   std::string corrLevel_;
00066   edm::InputTag jets_;
00068   std::map<std::string,TH1F*> hists_; 
00069 };
00070 
00071 
00072 PatJetAnalyzer::PatJetAnalyzer(const edm::ParameterSet& cfg):
00073   corrLevel_(cfg.getParameter<std::string>("corrLevel")),
00074   jets_(cfg.getParameter<edm::InputTag>("src"))
00075 {
00076   // register TFileService
00077   edm::Service<TFileService> fs;
00078 
00079   // jet multiplicity
00080   hists_["mult" ]=fs->make<TH1F>("mult" , "N_{Jet}"          ,   15,   0.,   15.);
00081   // jet pt (for all jets)
00082   hists_["pt"   ]=fs->make<TH1F>("pt"   , "p_{T}(Jet) [GeV]" ,   60,   0.,  300.);
00083   // jet eta (for all jets)
00084   hists_["eta"  ]=fs->make<TH1F>("eta"  , "#eta (Jet)"       ,   60,  -3.,    3.);
00085   // jet phi (for all jets)
00086   hists_["phi"  ]=fs->make<TH1F>("phi"  , "#phi (Jet)"       ,   60,  3.2,   3.2);
00087   // dijet mass (if available)
00088   hists_["mass" ]=fs->make<TH1F>("mass" , "M_{jj} [GeV]"     ,   50,   0.,  500.);
00089   // basic histograms for jet energy response
00090   for(unsigned int idx=0; idx<MAXBIN; ++idx){
00091     char buffer [10]; sprintf (buffer, "jes_%i", idx);
00092     char title  [50]; sprintf (title , "p_{T}^{rec}/p_{T}^{gen} [%i GeV - %i GeV]", (int)BINS[idx], (int)BINS[idx+1]);
00093     hists_[buffer]=fs->make<TH1F>(buffer, title,  100, 0., 2.);
00094   }  
00095 }
00096 
00097 void
00098 PatJetAnalyzer::analyze(const edm::Event& event, const edm::EventSetup& setup)
00099 {
00100   // recieve jet collection label
00101   edm::Handle<edm::View<pat::Jet> > jets;
00102   event.getByLabel(jets_,jets);
00103 
00104   // loop jets
00105   for(edm::View<pat::Jet>::const_iterator jet=jets->begin(); jet!=jets->end(); ++jet){
00106     // print jec factors
00107     // print(jet, jet-jets->begin());
00108 
00109     // fill basic kinematics
00110     fill( "pt" , jet->correctedJet(corrLevel_).pt());
00111     fill( "eta", jet->eta());
00112     fill( "phi", jet->phi());
00113     // basic plots for jet responds plot as a function of pt
00114     if( jet->genJet() ){
00115       double resp=jet->correctedJet(corrLevel_).pt()/jet->genJet()->pt();
00116       for(unsigned int idx=0; idx<MAXBIN; ++idx){
00117         if(BINS[idx]<=jet->genJet()->pt() && jet->genJet()->pt()<BINS[idx+1]){
00118           char buffer [10]; sprintf (buffer, "jes_%i", idx);
00119           fill( buffer, resp );
00120         }
00121       }
00122     }
00123   }
00124   // jet multiplicity
00125   fill( "mult" , jets->size());
00126   // invariant dijet mass for first two leading jets
00127   if(jets->size()>1){ fill( "mass", ((*jets)[0].p4()+(*jets)[1].p4()).mass());}
00128 }
00129 
00130 void
00131 PatJetAnalyzer::print(edm::View<pat::Jet>::const_iterator& jet, unsigned int idx)
00132 {
00133   //edm::LogInfo log("JEC");
00134   std::cout << "[" << idx << "] :: eta=" << std::setw(10) << jet->eta() << " phi=" << std::setw(10) << jet->phi() << " size: " << jet->availableJECLevels().size() << std::endl;
00135   for(unsigned int idx=0; idx<jet->availableJECLevels().size(); ++idx){
00136     pat::Jet correctedJet;
00137     if(jet->availableJECLevels()[idx].find("L5Flavor")!=std::string::npos|| 
00138        jet->availableJECLevels()[idx].find("L7Parton")!=std::string::npos ){
00139       correctedJet=jet->correctedJet(idx, pat::JetCorrFactors::UDS);
00140     }
00141     else{
00142       correctedJet=jet->correctedJet(idx, pat::JetCorrFactors::NONE );
00143     }
00144     std::cout << std::setw(10) << correctedJet.currentJECLevel() << " pt=" << std::setw(10) << correctedJet.pt() << std::endl;
00145   }
00146 }
00147 
00148 #include "FWCore/Framework/interface/MakerMacros.h"
00149 DEFINE_FWK_MODULE(PatJetAnalyzer);