CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
GenTrackMatcher.cc
Go to the documentation of this file.
1 
14 
15 namespace edm
16 {
17 class ParameterSet;
18 }
19 
21 {
22 public:
25 
26 private:
27  void produce( edm::Event& evt, const edm::EventSetup& es ) override;
31 };
32 
38 
39 using namespace edm;
40 using namespace std;
41 using namespace reco;
42 
44  tracer_(p),
45  tracks_(p.getUntrackedParameter<edm::InputTag>("trackProducer")),
46  genParticles_(p.getUntrackedParameter<edm::InputTag>("genParticles"))
47 {
48  produces<GenParticleMatch>();
49 }
50 
52 {
54  evt.getByLabel(tracks_, tracks);
55  Handle<vector<int> > barCodes;
56  evt.getByLabel(genParticles_, barCodes);
58  evt.getByLabel(genParticles_, genParticles);
59  auto_ptr<GenParticleMatch> match(new GenParticleMatch(GenParticleRefProd(genParticles)));
60  GenParticleMatch::Filler filler(*match);
61  size_t n = tracks->size();
62  vector<int> indices(n,-1);
63  tracer_.newEvent(evt, es);
64  for (size_t i = 0; i < n; ++ i )
65  {
66  RefToBase<Track> track(tracks, i);
67  if (tracer_.evaluate(track))
68  {
69  const HepMC::GenParticle * particle = tracer_.genParticle();
70  if (particle)
71  {
72  int barCode = particle->barcode();
73  vector<int>::const_iterator b = barCodes->begin(), e = barCodes->end(), f = find( b, e, barCode );
75  << "found matching particle with barcode" << *f
76  << " which has not been found in " << genParticles_;
77  indices[i] = *f;
78  }
79  }
80  }
81  filler.insert(tracks, indices.begin(), indices.end());
82  filler.fill();
83  evt.put(match);
84 }
85 
87 
89 
int i
Definition: DBlmapReader.cc:9
void newEvent(const edm::Event &, const edm::EventSetup &)
Pre-process event information (for accessing reconstruction information)
Definition: TrackHistory.cc:32
void produce(edm::Event &evt, const edm::EventSetup &es) override
const std::vector< reco::PFCandidatePtr > & tracks_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:52
edm::RefProd< GenParticleCollection > GenParticleRefProd
persistent reference to a GenParticle collection
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
edm::InputTag genParticles_
GenTrackMatcher(const edm::ParameterSet &)
constructor
double f[11][100]
This class traces the simulated and generated history of a given track.
Definition: TrackHistory.h:16
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:390
tuple tracks
Definition: testEve_cfg.py:39
double b
Definition: hdecay.h:120
edm::InputTag tracks_
TrackHistory tracer_
bool evaluate(TrackingParticleRef tpr)
Evaluate track history using a TrackingParticleRef.
Definition: TrackHistory.h:39
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:6
edm::Association< reco::GenParticleCollection > GenParticleMatch
const HepMC::GenParticle * genParticle() const
Returns a pointer to most primitive status 1 or 2 particle.
Definition: HistoryBase.h:89