File indexing completed on 2024-04-06 12:31:14
0001 #include "FWCore/Framework/interface/Frameworkfwd.h"
0002 #include "FWCore/Framework/interface/global/EDProducer.h"
0003 #include "FWCore/Framework/interface/Event.h"
0004 #include "FWCore/Utilities/interface/InputTag.h"
0005 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0006 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
0007
0008 namespace TopInitID {
0009 static constexpr int tID = 6;
0010 }
0011
0012 class TopInitSubset : public edm::global::EDProducer<> {
0013 public:
0014 explicit TopInitSubset(const edm::ParameterSet&);
0015
0016 void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
0017 void fillOutput(const reco::GenParticleCollection&, reco::GenParticleCollection&) const;
0018
0019 private:
0020 edm::EDGetTokenT<reco::GenParticleCollection> srcToken_;
0021 };
0022
0023 TopInitSubset::TopInitSubset(const edm::ParameterSet& cfg)
0024 : srcToken_(consumes<reco::GenParticleCollection>(cfg.getParameter<edm::InputTag>("src"))) {
0025 produces<reco::GenParticleCollection>();
0026 }
0027
0028 void TopInitSubset::produce(edm::StreamID, edm::Event& evt, const edm::EventSetup& setup) const {
0029 edm::Handle<reco::GenParticleCollection> src;
0030 evt.getByToken(srcToken_, src);
0031
0032 const reco::GenParticleRefProd ref = evt.getRefBeforePut<reco::GenParticleCollection>();
0033 auto sel = std::make_unique<reco::GenParticleCollection>();
0034
0035
0036 fillOutput(*src, *sel);
0037
0038 evt.put(std::move(sel));
0039 }
0040
0041 void TopInitSubset::fillOutput(const reco::GenParticleCollection& src, reco::GenParticleCollection& sel) const {
0042 for (auto const& t : src) {
0043 if (std::abs(t.pdgId()) == TopInitID::tID) {
0044 bool hasTopMother = false;
0045 for (unsigned idx = 0; idx < t.numberOfMothers(); ++idx)
0046 if (std::abs(t.mother(idx)->pdgId()) == TopInitID::tID)
0047 hasTopMother = true;
0048 if (hasTopMother)
0049 continue;
0050 for (unsigned idx = 0; idx < t.numberOfMothers(); ++idx) {
0051 sel.emplace_back(t.mother(idx)->threeCharge(),
0052 t.mother(idx)->p4(),
0053 t.mother(idx)->vertex(),
0054 t.mother(idx)->pdgId(),
0055 t.mother(idx)->status(),
0056 false);
0057 }
0058 break;
0059 }
0060 }
0061 }
0062
0063 #include "FWCore/Framework/interface/MakerMacros.h"
0064 DEFINE_FWK_MODULE(TopInitSubset);