File indexing completed on 2024-04-06 12:25:19
0001 #include <memory>
0002 #include <vector>
0003
0004 #include "FWCore/Framework/interface/Frameworkfwd.h"
0005 #include "FWCore/Framework/interface/global/EDProducer.h"
0006
0007 #include "FWCore/Framework/interface/Event.h"
0008 #include "FWCore/Framework/interface/MakerMacros.h"
0009
0010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0011 #include "FWCore/Utilities/interface/InputTag.h"
0012
0013 #include "DataFormats/Common/interface/View.h"
0014 #include "DataFormats/JetReco/interface/GenJetCollection.h"
0015 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
0016 #include "DataFormats/HepMCCandidate/interface/GenParticleFwd.h"
0017
0018 class HiSignalParticleProducer : public edm::global::EDProducer<> {
0019 public:
0020 explicit HiSignalParticleProducer(const edm::ParameterSet&);
0021 ~HiSignalParticleProducer() override = default;
0022
0023 static void fillDescriptions(edm::ConfigurationDescriptions&);
0024
0025 private:
0026 void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
0027
0028 edm::EDGetTokenT<edm::View<reco::GenParticle> > genParticleSrc_;
0029 };
0030
0031 HiSignalParticleProducer::HiSignalParticleProducer(const edm::ParameterSet& iConfig)
0032 : genParticleSrc_(consumes<edm::View<reco::GenParticle> >(iConfig.getParameter<edm::InputTag>("src"))) {
0033 std::string alias = (iConfig.getParameter<edm::InputTag>("src")).label();
0034 produces<reco::GenParticleCollection>().setBranchAlias(alias);
0035 }
0036
0037 void HiSignalParticleProducer::produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup&) const {
0038 auto signalGenParticles = std::make_unique<reco::GenParticleCollection>();
0039
0040 edm::Handle<edm::View<reco::GenParticle> > genParticles;
0041 iEvent.getByToken(genParticleSrc_, genParticles);
0042
0043 for (const reco::GenParticle& genParticle : *genParticles) {
0044 if (genParticle.collisionId() == 0) {
0045 signalGenParticles->push_back(genParticle);
0046 }
0047 }
0048
0049 iEvent.put(std::move(signalGenParticles));
0050 }
0051
0052 void HiSignalParticleProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0053 edm::ParameterSetDescription desc;
0054 desc.setComment("Selects genParticles from collision id = 0");
0055 desc.add<edm::InputTag>("src", edm::InputTag("genParticles"));
0056 descriptions.addWithDefaultLabel(desc);
0057 }
0058
0059 DEFINE_FWK_MODULE(HiSignalParticleProducer);