CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MCTrackMatcher.cc
Go to the documentation of this file.
1 
11 
12 namespace edm { class ParameterSet; }
13 
15  public:
18 
19  private:
20  void produce( edm::Event& evt, const edm::EventSetup& es ) override;
24 };
25 
34 using namespace edm;
35 using namespace std;
36 using namespace reco;
37 
39  associator_(p.getParameter<string>("associator")),
40  tracks_(p.getParameter<InputTag>("tracks")),
41  genParticles_( p.getParameter<InputTag>("genParticles")),
42  trackingParticles_( p.getParameter<InputTag>("trackingParticles")) {
43  produces<GenParticleMatch>();
44 }
45 
46 void MCTrackMatcher::produce(Event& evt, const EventSetup& es) {
48  es.get<TrackAssociatorRecord>().get(associator_,assoc);
49  const TrackAssociatorBase * associator = assoc.product();
51  evt.getByLabel(tracks_, tracks);
53  evt.getByLabel(trackingParticles_,trackingParticles);
54  Handle<vector<int> > barCodes;
55  evt.getByLabel(genParticles_,barCodes );
57  evt.getByLabel(genParticles_, genParticles );
58  RecoToSimCollection associations = associator->associateRecoToSim ( tracks, trackingParticles, & evt, &es );
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  for (size_t i = 0; i < n; ++ i ) {
64  RefToBase<Track> track(tracks, i);
65  RecoToSimCollection::const_iterator f = associations.find(track);
66  if ( f != associations.end() ) {
67  TrackingParticleRef tp = f->val.front().first;
68  TrackingParticle::genp_iterator j, b = tp->genParticle_begin(), e = tp->genParticle_end();
69  for( j = b; j != e; ++ j ) {
70  const reco::GenParticle * p = j->get();
71  if (p->status() == 1) {
72  indices[i] = j->key();
73  break;
74  }
75  }
76  }
77  }
78  filler.insert(tracks, indices.begin(), indices.end());
79  filler.fill();
80  evt.put(match);
81 }
82 
84 
86 
edm::Association< reco::GenParticleCollection > GenParticleMatch
int i
Definition: DBlmapReader.cc:9
const std::vector< reco::PFCandidatePtr > & tracks_
const_iterator end() const
last iterator over the map (read only)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
const_iterator find(const key_type &k) const
find element with specified reference key
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:52
edm::RefProd< GenParticleCollection > GenParticleRefProd
persistent reference to a GenParticle collection
std::string associator_
virtual int status() const GCC11_FINAL
status word
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
edm::InputTag genParticles_
int j
Definition: DBlmapReader.cc:9
double f[11][100]
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:390
virtual reco::RecoToSimCollection associateRecoToSim(edm::Handle< edm::View< reco::Track > > &tCH, edm::Handle< TrackingParticleCollection > &tPCH, const edm::Event *event, const edm::EventSetup *setup) const
compare reco to sim the handle of reco::Track and TrackingParticle collections
tuple tracks
Definition: testEve_cfg.py:39
edm::InputTag trackingParticles_
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:62
double b
Definition: hdecay.h:120
MCTrackMatcher(const edm::ParameterSet &)
constructor
void produce(edm::Event &evt, const edm::EventSetup &es) override
edm::InputTag tracks_
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:6