CMS 3D CMS Logo

MCTrackMatcher.cc
Go to the documentation of this file.
1 
13 
14 namespace edm { class ParameterSet; }
15 
16 using namespace edm;
17 using namespace std;
18 using namespace reco;
19 
21  public:
24 
25  private:
26  void produce( edm::StreamID, edm::Event& evt, const edm::EventSetup& es ) const override;
33 };
34 
40 
42  associator_(consumes<reco::TrackToTrackingParticleAssociator>(p.getParameter<string>("associator"))),
43  tracks_(consumes<edm::View<reco::Track>>(p.getParameter<InputTag>("tracks"))),
44  genParticles_(consumes<GenParticleCollection>(p.getParameter<InputTag>("genParticles"))),
45  genParticleInts_(consumes<std::vector<int>>(p.getParameter<InputTag>("genParticles"))),
46  trackingParticles_(consumes<TrackingParticleCollection>(p.getParameter<InputTag>("trackingParticles"))) {
47  produces<GenParticleMatch>();
48 }
49 
50 void MCTrackMatcher::produce(edm::StreamID, Event& evt, const EventSetup& es) const {
52  evt.getByToken(associator_,assoc);
55  evt.getByToken(tracks_, tracks);
57  evt.getByToken(trackingParticles_,trackingParticles);
58  Handle<vector<int> > barCodes;
59  evt.getByToken(genParticleInts_,barCodes );
61  evt.getByToken(genParticles_, genParticles );
62  RecoToSimCollection associations = associator->associateRecoToSim ( tracks, trackingParticles);
63  unique_ptr<GenParticleMatch> match(new GenParticleMatch(GenParticleRefProd(genParticles)));
65  size_t n = tracks->size();
66  vector<int> indices(n,-1);
67  for (size_t i = 0; i < n; ++ i ) {
68  RefToBase<Track> track(tracks, i);
69  RecoToSimCollection::const_iterator f = associations.find(track);
70  if ( f != associations.end() ) {
71  TrackingParticleRef tp = f->val.front().first;
72  TrackingParticle::genp_iterator j, b = tp->genParticle_begin(), e = tp->genParticle_end();
73  for( j = b; j != e; ++ j ) {
74  const reco::GenParticle * p = j->get();
75  if (p->status() == 1) {
76  indices[i] = j->key();
77  break;
78  }
79  }
80  }
81  }
82  filler.insert(tracks, indices.begin(), indices.end());
83  filler.fill();
84  evt.put(std::move(match));
85 }
86 
88 
90 
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
edm::Association< reco::GenParticleCollection > GenParticleMatch
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:127
const std::vector< reco::PFCandidatePtr > & tracks_
void produce(edm::StreamID, edm::Event &evt, const edm::EventSetup &es) const override
std::vector< TrackingParticle > TrackingParticleCollection
const_iterator end() const
last iterator over the map (read only)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:508
#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:53
edm::RefProd< GenParticleCollection > GenParticleRefProd
persistent reference to a GenParticle collection
edm::EDGetTokenT< reco::TrackToTrackingParticleAssociator > associator_
edm::EDGetTokenT< edm::View< reco::Track > > tracks_
double f[11][100]
edm::EDGetTokenT< std::vector< int > > genParticleInts_
T const * product() const
Definition: Handle.h:81
reco::RecoToSimCollection associateRecoToSim(const edm::Handle< edm::View< reco::Track > > &tCH, const edm::Handle< TrackingParticleCollection > &tPCH) const
compare reco to sim the handle of reco::Track and TrackingParticle collections
double b
Definition: hdecay.h:120
MCTrackMatcher(const edm::ParameterSet &)
constructor
fixed size matrix
HLT enums.
int status() const final
status word
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
edm::EDGetTokenT< GenParticleCollection > genParticles_
edm::EDGetTokenT< TrackingParticleCollection > trackingParticles_
def move(src, dest)
Definition: eostools.py:510