Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2021-02-14 13:34:55

0001 // Implementation of template class: JetCorExample
0002 // Description:  Example of simple EDAnalyzer correcting jets "on the fly".
0003 // Author: K. Kousouris
0004 // Date:  25 - August - 2008
0005 #include "RecoJets/JetAnalyzers/interface/JetCorExample.h"
0006 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
0007 #include "DataFormats/JetReco/interface/PFJetCollection.h"
0008 #include "DataFormats/JetReco/interface/GenJetCollection.h"
0009 #include "DataFormats/JetReco/interface/CaloJet.h"
0010 #include "DataFormats/JetReco/interface/PFJet.h"
0011 #include "DataFormats/JetReco/interface/GenJet.h"
0012 #include "FWCore/Framework/interface/Event.h"
0013 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0014 #include "JetMETCorrections/Objects/interface/JetCorrector.h"
0015 #include <TFile.h>
0016 #include <cmath>
0017 using namespace edm;
0018 using namespace reco;
0019 using namespace std;
0020 ////////////////////////////////////////////////////////////////////////////////////////
0021 template <class Jet>
0022 JetCorExample<Jet>::JetCorExample(edm::ParameterSet const& cfg) {
0023   JetAlgorithm = cfg.getParameter<std::string>("JetAlgorithm");
0024   HistoFileName = cfg.getParameter<std::string>("HistoFileName");
0025   JetCorrectionService = cfg.getParameter<std::string>("JetCorrectionService");
0026 }
0027 ////////////////////////////////////////////////////////////////////////////////////////
0028 template <class Jet>
0029 void JetCorExample<Jet>::beginJob() {
0030   TString hname;
0031   m_file = new TFile(HistoFileName.c_str(), "RECREATE");
0032   /////////// Booking histograms //////////////////////////
0033   hname = "JetPt";
0034   m_HistNames1D[hname] = new TH1F(hname, hname, 100, 0, 1000);
0035   hname = "CorJetPt";
0036   m_HistNames1D[hname] = new TH1F(hname, hname, 100, 0, 1000);
0037 }
0038 ////////////////////////////////////////////////////////////////////////////////////////
0039 template <class Jet>
0040 void JetCorExample<Jet>::analyze(edm::Event const& evt, edm::EventSetup const& iSetup) {
0041   /////////// Get the jet collection //////////////////////
0042   Handle<JetCollection> jets;
0043   evt.getByLabel(JetAlgorithm, jets);
0044   typename JetCollection::const_iterator i_jet;
0045   TString hname;
0046   const JetCorrector* corrector = JetCorrector::getJetCorrector(JetCorrectionService, iSetup);
0047   double scale;
0048   /////////// Loop over all jets and apply correction /////
0049   for (i_jet = jets->begin(); i_jet != jets->end(); i_jet++) {
0050     scale = corrector->correction(i_jet->p4());
0051     hname = "JetPt";
0052     FillHist1D(hname, i_jet->pt());
0053     hname = "CorJetPt";
0054     FillHist1D(hname, scale * i_jet->pt());
0055   }
0056 }
0057 ////////////////////////////////////////////////////////////////////////////////////////
0058 template <class Jet>
0059 void JetCorExample<Jet>::endJob() {
0060   /////////// Write Histograms in output ROOT file ////////
0061   if (m_file != nullptr) {
0062     m_file->cd();
0063     for (std::map<TString, TH1*>::iterator hid = m_HistNames1D.begin(); hid != m_HistNames1D.end(); hid++)
0064       hid->second->Write();
0065     delete m_file;
0066     m_file = nullptr;
0067   }
0068 }
0069 ////////////////////////////////////////////////////////////////////////////////////////
0070 template <class Jet>
0071 void JetCorExample<Jet>::FillHist1D(const TString& histName, const Double_t& value) {
0072   std::map<TString, TH1*>::iterator hid = m_HistNames1D.find(histName);
0073   if (hid == m_HistNames1D.end())
0074     std::cout << "%fillHist -- Could not find histogram with name: " << histName << std::endl;
0075   else
0076     hid->second->Fill(value);
0077 }
0078 /////////// Register Modules ////////
0079 #include "FWCore/Framework/interface/MakerMacros.h"
0080 /////////// Calo Jet Instance ////////
0081 typedef JetCorExample<CaloJet> CaloJetCorExample;
0082 DEFINE_FWK_MODULE(CaloJetCorExample);
0083 /////////// PF Jet Instance ////////
0084 typedef JetCorExample<PFJet> PFJetCorExample;
0085 DEFINE_FWK_MODULE(PFJetCorExample);