CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MuonLinksProducerForHLT.cc
Go to the documentation of this file.
1 
7 // system include files
8 #include <memory>
9 
10 // user include files
13 
16 
18 
21 
22 //#include <algorithm>
23 
25 {
26  produces<reco::MuonTrackLinksCollection>();
27  theLinkCollectionInInput = iConfig.getParameter<edm::InputTag>("LinkCollection");
28  theInclusiveTrackCollectionInInput = iConfig.getParameter<edm::InputTag>("InclusiveTrackerTrackCollection");
29  ptMin = iConfig.getParameter<double>("ptMin");
30  pMin = iConfig.getParameter<double>("pMin");
31  shareHitFraction = iConfig.getParameter<double>("shareHitFraction");
32 
33  linkToken_ = consumes<reco::MuonTrackLinksCollection>(theLinkCollectionInInput);
34  trackToken_ = consumes<reco::TrackCollection>(theInclusiveTrackCollectionInInput);
35 
36 }
37 
39 {
40 }
41 
43 {
44  std::auto_ptr<reco::MuonTrackLinksCollection> output(new reco::MuonTrackLinksCollection());
45 
47  iEvent.getByToken(linkToken_, links);
48 
50  iEvent.getByToken(trackToken_, incTracks);
51 
52  for(reco::MuonTrackLinksCollection::const_iterator link = links->begin();
53  link != links->end(); ++link){
54  bool found = false;
55  unsigned int trackIndex = 0;
56  unsigned int muonTrackHits = link->trackerTrack()->extra()->recHitsSize();
57  for(reco::TrackCollection::const_iterator track = incTracks->begin();
58  track != incTracks->end(); ++track, ++trackIndex){
59  if ( track->pt() < ptMin ) continue;
60  if ( track->p() < pMin ) continue;
61  //std::cout << "pt (muon/track) " << link->trackerTrack()->pt() << " " << track->pt() << std::endl;
62  unsigned trackHits = track->extra()->recHitsSize();
63  //std::cout << "hits (muon/track) " << muonTrackHits << " " << trackHits() << std::endl;
64  unsigned int smallestNumberOfHits = trackHits < muonTrackHits ? trackHits : muonTrackHits;
65  int numberOfCommonDetIds = 0;
66  for ( auto hit = track->extra()->recHitsBegin();
67  hit != track->extra()->recHitsEnd(); ++hit ) {
68  for ( auto mit = link->trackerTrack()->extra()->recHitsBegin();
69  mit != link->trackerTrack()->extra()->recHitsEnd(); ++mit ) {
70  if ( (*hit)->geographicalId() == (*mit)->geographicalId() &&
71  (*hit)->sharesInput((*mit),TrackingRecHit::some) ) {
72  numberOfCommonDetIds++;
73  break;
74  }
75  }
76  }
77  double fraction = (double)numberOfCommonDetIds/smallestNumberOfHits;
78  // std::cout << "Overlap/Smallest/fraction = " << numberOfCommonDetIds << " " << smallestNumberOfHits << " " << fraction << std::endl;
79  if( fraction > shareHitFraction ) {
80  output->push_back(reco::MuonTrackLinks(reco::TrackRef(incTracks,trackIndex),
81  link->standAloneTrack(),
82  link->globalTrack() ) );
83  found = true;
84  break;
85  }
86  }
87  if (!found)
88  output->push_back(reco::MuonTrackLinks(link->trackerTrack(),
89  link->standAloneTrack(),
90  link->globalTrack() ) );
91  }
92  iEvent.put( output );
93 }
T getParameter(std::string const &) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:446
std::vector< MuonTrackLinks > MuonTrackLinksCollection
collection of MuonTrackLinks
Definition: MuonFwd.h:22
int iEvent
Definition: GenABIO.cc:230
edm::EDGetTokenT< reco::TrackCollection > trackToken_
edm::EDGetTokenT< reco::MuonTrackLinksCollection > linkToken_
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
MuonLinksProducerForHLT(const edm::ParameterSet &)
edm::InputTag theInclusiveTrackCollectionInInput
virtual void produce(edm::Event &, const edm::EventSetup &)