CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10/src/RecoJets/JetAnalyzers/src/JetPlotsExample.cc

Go to the documentation of this file.
00001 // Implementation of template class: JetPlotsExample
00002 // Description:  Example of simple EDAnalyzer for jets.
00003 // Author: K. Kousouris
00004 // Date:  25 - August - 2008
00005 #include "RecoJets/JetAnalyzers/interface/JetPlotsExample.h"
00006 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00007 #include "DataFormats/JetReco/interface/PFJetCollection.h"
00008 #include "DataFormats/JetReco/interface/GenJetCollection.h"
00009 #include "DataFormats/JetReco/interface/CaloJet.h"
00010 #include "DataFormats/JetReco/interface/PFJet.h"
00011 #include "DataFormats/JetReco/interface/GenJet.h"
00012 #include "DataFormats/JetReco/interface/JPTJet.h"
00013 #include "DataFormats/JetReco/interface/JPTJetCollection.h"
00014 
00015 
00016 #include "FWCore/Framework/interface/Event.h"
00017 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00018 #include <TFile.h>
00019 #include <cmath>
00020 using namespace edm;
00021 using namespace reco;
00022 using namespace std;
00024 template<class Jet>
00025 JetPlotsExample<Jet>::JetPlotsExample(edm::ParameterSet const& cfg)
00026 {
00027   JetAlgorithm  = cfg.getParameter<std::string> ("JetAlgorithm"); 
00028   HistoFileName = cfg.getParameter<std::string> ("HistoFileName");
00029   NJets         = cfg.getParameter<int> ("NJets");
00030 }
00032 template<class Jet>
00033 void JetPlotsExample<Jet>::beginJob() 
00034 {
00035   TString hname;
00036   m_file = new TFile(HistoFileName.c_str(),"RECREATE"); 
00038   hname = "JetPt";
00039   m_HistNames1D[hname] = new TH1F(hname,hname,100,0,1000);
00040   hname = "JetEta";
00041   m_HistNames1D[hname] = new TH1F(hname,hname,120,-6,6);
00042   hname = "JetPhi";
00043   m_HistNames1D[hname] = new TH1F(hname,hname,100,-M_PI,M_PI);
00044   hname = "NumberOfJets";
00045   m_HistNames1D[hname] = new TH1F(hname,hname,100,0,100);
00046 }
00048 template<class Jet>
00049 void JetPlotsExample<Jet>::analyze(edm::Event const& evt, edm::EventSetup const& iSetup) 
00050 {
00052   Handle<JetCollection> jets;
00053   evt.getByLabel(JetAlgorithm,jets);
00054   typename JetCollection::const_iterator i_jet;
00055   int index = 0;
00056   TString hname; 
00058   hname = "NumberOfJets";
00059   FillHist1D(hname,jets->size()); 
00061   for(i_jet = jets->begin(); i_jet != jets->end() && index < NJets; ++i_jet) 
00062     {
00063       hname = "JetPt";
00064       FillHist1D(hname,i_jet->pt());   
00065       hname = "JetEta";
00066       FillHist1D(hname,i_jet->eta());
00067       hname = "JetPhi";
00068       FillHist1D(hname,i_jet->phi());
00069       index++;
00070     }
00071 }
00073 template<class Jet>
00074 void JetPlotsExample<Jet>::endJob() 
00075 {
00077   if (m_file !=0) 
00078     {
00079       m_file->cd();
00080       for (std::map<TString, TH1*>::iterator hid = m_HistNames1D.begin(); hid != m_HistNames1D.end(); hid++)
00081         hid->second->Write();
00082       delete m_file;
00083       m_file = 0;      
00084     }
00085 }
00087 template<class Jet>
00088 void JetPlotsExample<Jet>::FillHist1D(const TString& histName,const Double_t& value) 
00089 {
00090   std::map<TString, TH1*>::iterator hid=m_HistNames1D.find(histName);
00091   if (hid==m_HistNames1D.end())
00092     std::cout << "%fillHist -- Could not find histogram with name: " << histName << std::endl;
00093   else
00094     hid->second->Fill(value);
00095 }
00097 #include "FWCore/Framework/interface/MakerMacros.h"
00099 typedef JetPlotsExample<CaloJet> CaloJetPlotsExample;
00100 DEFINE_FWK_MODULE(CaloJetPlotsExample);
00102 typedef JetPlotsExample<GenJet> GenJetPlotsExample;
00103 DEFINE_FWK_MODULE(GenJetPlotsExample);
00105 typedef JetPlotsExample<PFJet> PFJetPlotsExample;
00106 DEFINE_FWK_MODULE(PFJetPlotsExample);
00108 typedef JetPlotsExample<JPTJet> JPTJetPlotsExample;
00109 DEFINE_FWK_MODULE(JPTJetPlotsExample);