CMS 3D CMS Logo

TrackingParticleBHadronRefSelector.cc
Go to the documentation of this file.
7 
11 
13 
14 #include "HepPDT/ParticleID.hh"
15 
17 public:
19 
20  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
21 
22  void produce(edm::Event& iEvent, const edm::EventSetup& iSetup) override;
23 
24 private:
26 
28 };
29 
30 
32  tpToken_(consumes<TrackingParticleCollection>(iConfig.getParameter<edm::InputTag>("src")))
33 {
34  tracer_.depth(-2); // as in SimTracker/TrackHistory/src/TrackClassifier.cc
35 
36  produces<TrackingParticleRefVector>();
37 }
38 
41  desc.add<edm::InputTag>("src", edm::InputTag("mix", "MergedTrackTruth"));
42  descriptions.add("trackingParticleBHadronRefSelector", desc);
43 }
44 
47  iEvent.getByToken(tpToken_, h_tps);
48 
49  auto ret = std::make_unique<TrackingParticleRefVector>();
50 
51  // Logic is similar to SimTracker/TrackHistory
52  for(size_t i=0, end=h_tps->size(); i<end; ++i) {
53  auto tpRef = TrackingParticleRef(h_tps, i);
54  if(tracer_.evaluate(tpRef)) { // ignore TP if history can not be traced
55  // following is from TrackClassifier::processesAtGenerator()
56  HistoryBase::RecoGenParticleTrail const & recoGenParticleTrail = tracer_.recoGenParticleTrail();
57  for(const auto& particle: recoGenParticleTrail) {
58  HepPDT::ParticleID particleID(particle->pdgId());
59  if(particleID.hasBottom()) {
60  ret->push_back(tpRef);
61  break;
62  }
63  }
64  }
65  }
66 
67  iEvent.put(std::move(ret));
68 }
69 
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:127
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
std::vector< TrackingParticle > TrackingParticleCollection
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:508
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
std::vector< const reco::GenParticle * > RecoGenParticleTrail
reco::GenParticle trail type.
Definition: HistoryBase.h:21
int iEvent
Definition: GenABIO.cc:230
edm::EDGetTokenT< TrackingParticleCollection > tpToken_
#define end
Definition: vmac.h:37
ParameterDescriptionBase * add(U const &iLabel, T const &value)
TrackingParticleBHadronRefSelector(const edm::ParameterSet &iConfig)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Base class to all the history types.
Definition: HistoryBase.h:12
bool evaluate(TrackingParticleRef tpr)
Evaluate track history using a TrackingParticleRef.
Definition: HistoryBase.h:122
HLT enums.
void depth(int d)
Set the depth of the history.
Definition: HistoryBase.h:53
RecoGenParticleTrail const & recoGenParticleTrail() const
Return all reco::GenParticle in the history.
Definition: HistoryBase.h:83
edm::Ref< TrackingParticleCollection > TrackingParticleRef
def move(src, dest)
Definition: eostools.py:510