CMS 3D CMS Logo

MuonLinksProducerForHLT.cc
Go to the documentation of this file.
1 
6 // system include files
7 #include <memory>
8 
9 // user include files
11 
14 
16 
19 
20 //#include <algorithm>
21 
23  produces<reco::MuonTrackLinksCollection>();
24  theLinkCollectionInInput = iConfig.getParameter<edm::InputTag>("LinkCollection");
25  theInclusiveTrackCollectionInInput = iConfig.getParameter<edm::InputTag>("InclusiveTrackerTrackCollection");
26  ptMin = iConfig.getParameter<double>("ptMin");
27  pMin = iConfig.getParameter<double>("pMin");
28  shareHitFraction = iConfig.getParameter<double>("shareHitFraction");
29 
30  linkToken_ = consumes<reco::MuonTrackLinksCollection>(theLinkCollectionInInput);
31  trackToken_ = consumes<reco::TrackCollection>(theInclusiveTrackCollectionInInput);
32 }
33 
35 
37  auto output = std::make_unique<reco::MuonTrackLinksCollection>();
38 
40  iEvent.getByToken(linkToken_, links);
41 
43  iEvent.getByToken(trackToken_, incTracks);
44 
45  for (reco::MuonTrackLinksCollection::const_iterator link = links->begin(); link != links->end(); ++link) {
46  bool found = false;
47  unsigned int trackIndex = 0;
48  unsigned int muonTrackHits = link->trackerTrack()->extra()->recHitsSize();
49  for (reco::TrackCollection::const_iterator track = incTracks->begin(); track != incTracks->end();
50  ++track, ++trackIndex) {
51  if (track->pt() < ptMin)
52  continue;
53  if (track->p() < pMin)
54  continue;
55  //std::cout << "pt (muon/track) " << link->trackerTrack()->pt() << " " << track->pt() << std::endl;
56  unsigned trackHits = track->extra()->recHitsSize();
57  //std::cout << "hits (muon/track) " << muonTrackHits << " " << trackHits() << std::endl;
58  unsigned int smallestNumberOfHits = trackHits < muonTrackHits ? trackHits : muonTrackHits;
59  int numberOfCommonDetIds = 0;
60  for (auto hit = track->extra()->recHitsBegin(); hit != track->extra()->recHitsEnd(); ++hit) {
61  for (auto mit = link->trackerTrack()->extra()->recHitsBegin();
62  mit != link->trackerTrack()->extra()->recHitsEnd();
63  ++mit) {
64  if ((*hit)->geographicalId() == (*mit)->geographicalId() &&
65  (*hit)->sharesInput((*mit), TrackingRecHit::some)) {
66  numberOfCommonDetIds++;
67  break;
68  }
69  }
70  }
71  double fraction = (double)numberOfCommonDetIds / smallestNumberOfHits;
72  // std::cout << "Overlap/Smallest/fraction = " << numberOfCommonDetIds << " " << smallestNumberOfHits << " " << fraction << std::endl;
73  if (fraction > shareHitFraction) {
74  output->push_back(
75  reco::MuonTrackLinks(reco::TrackRef(incTracks, trackIndex), link->standAloneTrack(), link->globalTrack()));
76  found = true;
77  break;
78  }
79  }
80  if (!found)
81  output->push_back(reco::MuonTrackLinks(link->trackerTrack(), link->standAloneTrack(), link->globalTrack()));
82  }
83  iEvent.put(std::move(output));
84 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
int iEvent
Definition: GenABIO.cc:224
edm::EDGetTokenT< reco::TrackCollection > trackToken_
edm::EDGetTokenT< reco::MuonTrackLinksCollection > linkToken_
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
MuonLinksProducerForHLT(const edm::ParameterSet &)
edm::InputTag theInclusiveTrackCollectionInInput
def move(src, dest)
Definition: eostools.py:511