File indexing completed on 2024-04-06 12:31:06
0001 #include "FWCore/Framework/interface/Event.h"
0002 #include "FWCore/Framework/interface/MakerMacros.h"
0003 #include "FWCore/Framework/interface/stream/EDProducer.h"
0004 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0005 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0006 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0007
0008 #include "DataFormats/Common/interface/EDProductfwd.h"
0009 #include "SimDataFormats/TrackingAnalysis/interface/TrackingParticle.h"
0010 #include "SimDataFormats/TrackingAnalysis/interface/TrackingParticleFwd.h"
0011
0012 #include "SimTracker/TrackHistory/interface/HistoryBase.h"
0013
0014 #include "HepPDT/ParticleID.hh"
0015
0016 class TrackingParticleBHadronRefSelector : public edm::stream::EDProducer<> {
0017 public:
0018 TrackingParticleBHadronRefSelector(const edm::ParameterSet &iConfig);
0019
0020 static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
0021
0022 void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override;
0023
0024 private:
0025 edm::EDGetTokenT<TrackingParticleCollection> tpToken_;
0026
0027 HistoryBase tracer_;
0028 };
0029
0030 TrackingParticleBHadronRefSelector::TrackingParticleBHadronRefSelector(const edm::ParameterSet &iConfig)
0031 : tpToken_(consumes<TrackingParticleCollection>(iConfig.getParameter<edm::InputTag>("src"))) {
0032 tracer_.depth(-2);
0033
0034 produces<TrackingParticleRefVector>();
0035 }
0036
0037 void TrackingParticleBHadronRefSelector::fillDescriptions(edm::ConfigurationDescriptions &descriptions) {
0038 edm::ParameterSetDescription desc;
0039 desc.add<edm::InputTag>("src", edm::InputTag("mix", "MergedTrackTruth"));
0040 descriptions.add("trackingParticleBHadronRefSelectorDefault", desc);
0041 }
0042
0043 void TrackingParticleBHadronRefSelector::produce(edm::Event &iEvent, const edm::EventSetup &iSetup) {
0044 edm::Handle<TrackingParticleCollection> h_tps;
0045 iEvent.getByToken(tpToken_, h_tps);
0046
0047 auto ret = std::make_unique<TrackingParticleRefVector>();
0048
0049
0050 for (size_t i = 0, end = h_tps->size(); i < end; ++i) {
0051 auto tpRef = TrackingParticleRef(h_tps, i);
0052 if (tracer_.evaluate(tpRef)) {
0053
0054 HistoryBase::RecoGenParticleTrail const &recoGenParticleTrail = tracer_.recoGenParticleTrail();
0055 for (const auto &particle : recoGenParticleTrail) {
0056 HepPDT::ParticleID particleID(particle->pdgId());
0057 if (particleID.hasBottom()) {
0058 ret->push_back(tpRef);
0059 break;
0060 }
0061 }
0062 }
0063 }
0064
0065 iEvent.put(std::move(ret));
0066 }
0067
0068 DEFINE_FWK_MODULE(TrackingParticleBHadronRefSelector);