File indexing completed on 2024-04-06 12:20:44
0001 #include "FWCore/Framework/interface/Frameworkfwd.h"
0002 #include "FWCore/Framework/interface/stream/EDProducer.h"
0003 #include "FWCore/Framework/interface/Event.h"
0004 #include "FWCore/Framework/interface/MakerMacros.h"
0005 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0006
0007 #include "DataFormats/L1THGCal/interface/HGCalMulticluster.h"
0008 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
0009 #include "DataFormats/Math/interface/deltaR.h"
0010
0011 namespace l1t {
0012 class HGC3DClusterGenMatchSelector : public edm::stream::EDProducer<> {
0013 public:
0014 explicit HGC3DClusterGenMatchSelector(const edm::ParameterSet &);
0015 ~HGC3DClusterGenMatchSelector() override {}
0016
0017 private:
0018 edm::EDGetTokenT<l1t::HGCalMulticlusterBxCollection> src_;
0019 edm::EDGetToken genParticleSrc_;
0020 double dR_;
0021 void produce(edm::Event &, const edm::EventSetup &) override;
0022
0023 };
0024 }
0025
0026 l1t::HGC3DClusterGenMatchSelector::HGC3DClusterGenMatchSelector(const edm::ParameterSet &iConfig)
0027 : src_(consumes<l1t::HGCalMulticlusterBxCollection>(iConfig.getParameter<edm::InputTag>("src"))),
0028 genParticleSrc_(consumes<reco::GenParticleCollection>(iConfig.getParameter<edm::InputTag>("genSrc"))),
0029 dR_(iConfig.getParameter<double>("dR")) {
0030 produces<l1t::HGCalMulticlusterBxCollection>();
0031 }
0032
0033 void l1t::HGC3DClusterGenMatchSelector::produce(edm::Event &iEvent, const edm::EventSetup &) {
0034 auto out = std::make_unique<l1t::HGCalMulticlusterBxCollection>();
0035
0036 edm::Handle<l1t::HGCalMulticlusterBxCollection> multiclusters;
0037 iEvent.getByToken(src_, multiclusters);
0038
0039 edm::Handle<reco::GenParticleCollection> genParticles;
0040 iEvent.getByToken(genParticleSrc_, genParticles);
0041
0042 for (int bx = multiclusters->getFirstBX(); bx <= multiclusters->getLastBX(); ++bx) {
0043 for (auto it = multiclusters->begin(bx), ed = multiclusters->end(bx); it != ed; ++it) {
0044 const auto &multicluster = *it;
0045 for (const auto &particle : *genParticles) {
0046 if (particle.status() != 1)
0047 continue;
0048 if (deltaR(multicluster, particle) < dR_) {
0049 out->push_back(bx, multicluster);
0050 break;
0051 }
0052 }
0053 }
0054 }
0055
0056 iEvent.put(std::move(out));
0057 }
0058 using l1t::HGC3DClusterGenMatchSelector;
0059 DEFINE_FWK_MODULE(HGC3DClusterGenMatchSelector);