File indexing completed on 2024-04-06 12:23:47
0001
0002 #include "FWCore/Framework/interface/stream/EDProducer.h"
0003 #include "DataFormats/Common/interface/RefVector.h"
0004
0005 #include "CommonTools/UtilAlgos/interface/StringCutObjectSelector.h"
0006 #include "CommonTools/UtilAlgos/interface/SingleObjectSelector.h"
0007 #include "CommonTools/UtilAlgos/interface/ObjectSelector.h"
0008 #include "CommonTools/UtilAlgos/interface/SingleElementCollectionSelector.h"
0009
0010 #include "DataFormats/JetReco/interface/CaloJet.h"
0011 #include "DataFormats/Candidate/interface/Candidate.h"
0012
0013 #include <vector>
0014
0015 class CaloJetSlimmer : public edm::stream::EDProducer<> {
0016 public:
0017 CaloJetSlimmer(edm::ParameterSet const& params)
0018 : srcToken_(consumes<edm::View<reco::CaloJet> >(params.getParameter<edm::InputTag>("src"))),
0019 cut_(params.getParameter<std::string>("cut")),
0020 selector_(cut_) {
0021 produces<reco::CaloJetCollection>();
0022 }
0023
0024 ~CaloJetSlimmer() override {}
0025
0026 void produce(edm::Event& iEvent, const edm::EventSetup& iSetup) override {
0027 auto caloJets = std::make_unique<reco::CaloJetCollection>();
0028
0029 edm::Handle<edm::View<reco::CaloJet> > h_jets;
0030 iEvent.getByToken(srcToken_, h_jets);
0031
0032 for (auto const& ijet : *h_jets) {
0033 if (selector_(ijet)) {
0034 reco::CaloJet::Specific tmp_specific;
0035 const reco::Candidate::Point orivtx(0, 0, 0);
0036 tmp_specific.mEnergyFractionEm = ijet.emEnergyFraction();
0037 tmp_specific.mEnergyFractionHadronic = ijet.energyFractionHadronic();
0038 reco::CaloJet newCaloJet(ijet.p4(), orivtx, tmp_specific);
0039 float jetArea = ijet.jetArea();
0040 newCaloJet.setJetArea(jetArea);
0041 caloJets->push_back(newCaloJet);
0042 }
0043 }
0044 iEvent.put(std::move(caloJets), "");
0045 }
0046
0047 protected:
0048 const edm::EDGetTokenT<edm::View<reco::CaloJet> > srcToken_;
0049 const std::string cut_;
0050 const StringCutObjectSelector<reco::Jet> selector_;
0051 };
0052
0053 #include "FWCore/Framework/interface/MakerMacros.h"
0054 DEFINE_FWK_MODULE(CaloJetSlimmer);