Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 10:53:46

0001 // Package:    SinglePhotonJetPlusHOFilter
0002 // Class:      SinglePhotonJetPlusHOFilter
0003 //
0004 /*
0005  Description: [one line class summary]
0006 Skimming of SinglePhoton data set for the study of HO absolute weight calculation
0007 * Skimming Efficiency : ~ 2 %
0008  Implementation:
0009      [Notes on implementation]
0010      For Secondary Datasets (SD)
0011 */
0012 //
0013 // Original Author:  Gobinda Majumder & Suman Chatterjee
0014 //         Created:  Fri July 29 14:52:17 IST 2016
0015 // $Id$
0016 //
0017 //
0018 
0019 // system include files
0020 #include <memory>
0021 
0022 // class declaration
0023 //
0024 
0025 #include "FWCore/Framework/interface/global/EDFilter.h"
0026 
0027 #include "FWCore/Framework/interface/Event.h"
0028 #include "FWCore/Framework/interface/MakerMacros.h"
0029 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0030 
0031 #include "DataFormats/JetReco/interface/PFJet.h"
0032 #include "DataFormats/ParticleFlowReco/interface/PFClusterFwd.h"
0033 #include "DataFormats/ParticleFlowReco/interface/PFCluster.h"
0034 
0035 #include "DataFormats/EgammaCandidates/interface/Photon.h"
0036 #include "DataFormats/EgammaCandidates/interface/PhotonFwd.h"
0037 #include "DataFormats/Math/interface/deltaR.h"
0038 
0039 class SinglePhotonJetPlusHOFilter : public edm::global::EDFilter<> {
0040 public:
0041   explicit SinglePhotonJetPlusHOFilter(const edm::ParameterSet&);
0042 
0043   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0044 
0045 private:
0046   bool filter(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
0047 
0048   // ----------member data ---------------------------
0049   const double jtptthr_;
0050   const double jtetath_;
0051   const double hothres_;
0052   const double pho_Ptcut_;
0053 
0054   const edm::EDGetTokenT<reco::PFJetCollection> tok_PFJets_;
0055   const edm::EDGetTokenT<reco::PFClusterCollection> tok_hoht_;
0056   const edm::EDGetTokenT<edm::View<reco::Photon> > tok_photons_;
0057 };
0058 
0059 SinglePhotonJetPlusHOFilter::SinglePhotonJetPlusHOFilter(const edm::ParameterSet& iConfig)
0060     : jtptthr_{iConfig.getParameter<double>("Ptcut")},
0061       jtetath_{iConfig.getParameter<double>("Etacut")},
0062       hothres_{iConfig.getParameter<double>("HOcut")},
0063       pho_Ptcut_{iConfig.getParameter<double>("Pho_Ptcut")},
0064       tok_PFJets_{consumes<reco::PFJetCollection>(iConfig.getParameter<edm::InputTag>("PFJets"))},
0065       tok_hoht_{consumes<reco::PFClusterCollection>(iConfig.getParameter<edm::InputTag>("particleFlowClusterHO"))},
0066       tok_photons_{consumes<edm::View<reco::Photon> >(iConfig.getParameter<edm::InputTag>("Photons"))} {}
0067 
0068 // ------------ method called on each new Event  ------------
0069 bool SinglePhotonJetPlusHOFilter::filter(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const {
0070   using namespace std;
0071   using namespace edm;
0072   using namespace reco;
0073 
0074   bool passed = false;
0075   vector<pair<double, double> > jetdirection;
0076   vector<double> jetspt;
0077   if (auto PFJets = iEvent.getHandle(tok_PFJets_)) {
0078     for (const auto& jet : *PFJets) {
0079       if ((jet.pt() < jtptthr_) || (abs(jet.eta()) > jtetath_))
0080         continue;
0081 
0082       jetdirection.emplace_back(jet.eta(), jet.phi());
0083       jetspt.push_back(jet.pt());
0084       passed = true;
0085     }
0086   }
0087   if (!passed)
0088     return passed;
0089 
0090   bool pho_passed = false;
0091   vector<pair<double, double> > phodirection;
0092   vector<double> phopT;
0093   if (auto photons = iEvent.getHandle(tok_photons_)) {
0094     for (const auto& gamma1 : *photons) {
0095       if (gamma1.pt() < pho_Ptcut_)
0096         continue;
0097       phodirection.emplace_back(gamma1.eta(), gamma1.phi());
0098       phopT.push_back(gamma1.pt());
0099 
0100       pho_passed = true;
0101     }
0102   }
0103 
0104   if (!pho_passed)
0105     return false;
0106 
0107   bool isJetDir = false;
0108   if (auto hhoht = iEvent.getHandle(tok_hoht_)) {
0109     const auto& hoht = *hhoht;
0110     if (!hoht.empty()) {
0111       for (const auto& jet : jetdirection) {
0112         bool matched = false;
0113         for (const auto& ph : phodirection) {
0114           if (abs(deltaPhi(ph.second, jet.second)) > 2.0) {
0115             matched = true;
0116             break;
0117           }
0118         }
0119         if (matched) {
0120           for (const auto& ij : hoht) {
0121             double hoenr = ij.energy();
0122             if (hoenr < hothres_)
0123               continue;
0124 
0125             const math::XYZPoint& cluster_pos = ij.position();
0126 
0127             double hoeta = cluster_pos.eta();
0128             double hophi = cluster_pos.phi();
0129 
0130             double delta = deltaR2(jet.first, jet.second, hoeta, hophi);
0131             if (delta < 0.5) {
0132               isJetDir = true;
0133               break;
0134             }
0135           }
0136         }
0137         if (isJetDir) {
0138           break;
0139         }
0140       }
0141     }
0142   }
0143 
0144   return isJetDir;
0145 }
0146 
0147 // ------------ method fills 'descriptions' with the allowed parameters for the module  ------------
0148 void SinglePhotonJetPlusHOFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0149   //The following says we do not know what parameters are allowed so do no validation
0150   // Please change this to state exactly what you do use, even if it is no parameters
0151   edm::ParameterSetDescription desc;
0152   desc.add<double>("Ptcut", 90.0);
0153   desc.add<double>("Etacut", 1.5);
0154   desc.add<double>("HOcut", 5);
0155   desc.add<double>("Pho_Ptcut", 120);
0156   desc.add<edm::InputTag>("PFJets");
0157   desc.add<edm::InputTag>("particleFlowClusterHO");
0158   desc.add<edm::InputTag>("Photons");
0159   descriptions.addWithDefaultLabel(desc);
0160 }
0161 
0162 //define this as a plug-in
0163 DEFINE_FWK_MODULE(SinglePhotonJetPlusHOFilter);