Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2022-05-18 03:27:24

0001 // -*- C++ -*-
0002 //
0003 // Package:    PhysicsTools/NanoAOD
0004 // Class:      NanoAODBaseCrossCleaner
0005 //
0006 /**\class NanoAODBaseCrossCleaner NanoAODBaseCrossCleaner.cc PhysicsTools/NanoAODBaseCrossCleaner/plugins/NanoAODBaseCrossCleaner.cc
0007 
0008  Description: [one line class summary]
0009 
0010  Implementation:
0011      [Notes on implementation]
0012 */
0013 //
0014 // Original Author:  Andrea Rizzi
0015 //         Created:  Mon, 28 Aug 2017 09:26:39 GMT
0016 //
0017 //
0018 
0019 #include "PhysicsTools/NanoAOD/plugins/NanoAODBaseCrossCleaner.h"
0020 #include "DataFormats/NanoAOD/interface/FlatTable.h"
0021 
0022 //
0023 // constructors and destructor
0024 //
0025 NanoAODBaseCrossCleaner::NanoAODBaseCrossCleaner(const edm::ParameterSet& params)
0026     : name_(params.getParameter<std::string>("name")),
0027       doc_(params.getParameter<std::string>("doc")),
0028       jets_(consumes<edm::View<pat::Jet>>(params.getParameter<edm::InputTag>("jets"))),
0029       muons_(consumes<edm::View<pat::Muon>>(params.getParameter<edm::InputTag>("muons"))),
0030       electrons_(consumes<edm::View<pat::Electron>>(params.getParameter<edm::InputTag>("electrons"))),
0031       lowPtElectronsTag_(params.getParameter<edm::InputTag>("lowPtElectrons")),
0032       lowPtElectrons_(mayConsume<edm::View<pat::Electron>>(lowPtElectronsTag_)),
0033       taus_(consumes<edm::View<pat::Tau>>(params.getParameter<edm::InputTag>("taus"))),
0034       photons_(consumes<edm::View<pat::Photon>>(params.getParameter<edm::InputTag>("photons"))),
0035       jetSel_(params.getParameter<std::string>("jetSel")),
0036       muonSel_(params.getParameter<std::string>("muonSel")),
0037       electronSel_(params.getParameter<std::string>("electronSel")),
0038       lowPtElectronSel_(params.getParameter<std::string>("lowPtElectronSel")),
0039       tauSel_(params.getParameter<std::string>("tauSel")),
0040       photonSel_(params.getParameter<std::string>("photonSel")),
0041       jetName_(params.getParameter<std::string>("jetName")),
0042       muonName_(params.getParameter<std::string>("muonName")),
0043       electronName_(params.getParameter<std::string>("electronName")),
0044       lowPtElectronName_(params.getParameter<std::string>("lowPtElectronName")),
0045       tauName_(params.getParameter<std::string>("tauName")),
0046       photonName_(params.getParameter<std::string>("photonName"))
0047 
0048 {
0049   produces<nanoaod::FlatTable>("jets");
0050   produces<nanoaod::FlatTable>("muons");
0051   produces<nanoaod::FlatTable>("electrons");
0052   if (!lowPtElectronsTag_.label().empty())
0053     produces<nanoaod::FlatTable>("lowPtElectrons");
0054   produces<nanoaod::FlatTable>("taus");
0055   produces<nanoaod::FlatTable>("photons");
0056 }
0057 
0058 NanoAODBaseCrossCleaner::~NanoAODBaseCrossCleaner() {
0059   // do anything here that needs to be done at destruction time
0060   // (e.g. close files, deallocate resources etc.)
0061 }
0062 
0063 //
0064 // member functions
0065 //
0066 
0067 // ------------ method called to produce the data  ------------
0068 
0069 void NanoAODBaseCrossCleaner::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0070   using namespace edm;
0071   edm::Handle<edm::View<pat::Jet>> jetsIn;
0072   iEvent.getByToken(jets_, jetsIn);
0073   std::vector<uint8_t> jets;
0074   for (const auto& j : *jetsIn) {
0075     jets.push_back(jetSel_(j));
0076   }
0077   auto jetsTable = std::make_unique<nanoaod::FlatTable>(jetsIn->size(), jetName_, false, true);
0078 
0079   edm::Handle<edm::View<pat::Muon>> muonsIn;
0080   iEvent.getByToken(muons_, muonsIn);
0081   std::vector<uint8_t> muons;
0082   for (const auto& m : *muonsIn) {
0083     muons.push_back(muonSel_(m));
0084   }
0085   auto muonsTable = std::make_unique<nanoaod::FlatTable>(muonsIn->size(), muonName_, false, true);
0086 
0087   edm::Handle<edm::View<pat::Electron>> electronsIn;
0088   iEvent.getByToken(electrons_, electronsIn);
0089   std::vector<uint8_t> eles;
0090   for (const auto& e : *electronsIn) {
0091     eles.push_back(electronSel_(e));
0092   }
0093   auto electronsTable = std::make_unique<nanoaod::FlatTable>(electronsIn->size(), electronName_, false, true);
0094 
0095   edm::Handle<edm::View<pat::Electron>> lowPtElectronsIn;
0096   std::vector<uint8_t> lowPtEles;
0097   if (!lowPtElectronsTag_.label().empty()) {
0098     iEvent.getByToken(lowPtElectrons_, lowPtElectronsIn);
0099     for (const auto& e : *lowPtElectronsIn) {
0100       lowPtEles.push_back(lowPtElectronSel_(e));
0101     }
0102   }
0103 
0104   edm::Handle<edm::View<pat::Tau>> tausIn;
0105   iEvent.getByToken(taus_, tausIn);
0106   std::vector<uint8_t> taus;
0107   for (const auto& t : *tausIn) {
0108     taus.push_back(tauSel_(t));
0109   }
0110   auto tausTable = std::make_unique<nanoaod::FlatTable>(tausIn->size(), tauName_, false, true);
0111 
0112   edm::Handle<edm::View<pat::Photon>> photonsIn;
0113   iEvent.getByToken(photons_, photonsIn);
0114   std::vector<uint8_t> photons;
0115   for (const auto& p : *photonsIn) {
0116     photons.push_back(photonSel_(p));
0117   }
0118   auto photonsTable = std::make_unique<nanoaod::FlatTable>(photonsIn->size(), photonName_, false, true);
0119 
0120   objectSelection(*jetsIn, *muonsIn, *electronsIn, *tausIn, *photonsIn, jets, muons, eles, taus, photons);
0121 
0122   muonsTable->addColumn<uint8_t>(name_, muons, doc_);
0123   jetsTable->addColumn<uint8_t>(name_, jets, doc_);
0124   electronsTable->addColumn<uint8_t>(name_, eles, doc_);
0125   tausTable->addColumn<uint8_t>(name_, taus, doc_);
0126   photonsTable->addColumn<uint8_t>(name_, photons, doc_);
0127 
0128   iEvent.put(std::move(jetsTable), "jets");
0129   iEvent.put(std::move(muonsTable), "muons");
0130   iEvent.put(std::move(electronsTable), "electrons");
0131   iEvent.put(std::move(tausTable), "taus");
0132   iEvent.put(std::move(photonsTable), "photons");
0133 
0134   if (!lowPtElectronsTag_.label().empty()) {
0135     auto lowPtElectronsTable =
0136         std::make_unique<nanoaod::FlatTable>(lowPtElectronsIn->size(), lowPtElectronName_, false, true);
0137     lowPtElectronsTable->addColumn<uint8_t>(name_, lowPtEles, doc_);
0138     iEvent.put(std::move(lowPtElectronsTable), "lowPtElectrons");
0139   }
0140 }
0141 
0142 // ------------ method called once each stream before processing any runs, lumis or events  ------------
0143 void NanoAODBaseCrossCleaner::beginStream(edm::StreamID) {}
0144 
0145 // ------------ method called once each stream after processing all runs, lumis and events  ------------
0146 void NanoAODBaseCrossCleaner::endStream() {}
0147 
0148 // ------------ method fills 'descriptions' with the allowed parameters for the module  ------------
0149 void NanoAODBaseCrossCleaner::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0150   //The following says we do not know what parameters are allowed so do no validation
0151   // Please change this to state exactly what you do use, even if it is no parameters
0152   edm::ParameterSetDescription desc;
0153   desc.setUnknown();
0154   descriptions.addDefault(desc);
0155 }
0156 
0157 //define this as a plug-in
0158 DEFINE_FWK_MODULE(NanoAODBaseCrossCleaner);