CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CosmicMuonLinksProducer.cc
Go to the documentation of this file.
2 
9 // system include files
10 #include <memory>
11 
12 // user include files
14 
17 
19 
23 
27 
28 using namespace edm;
29 using namespace std;
30 using namespace reco;
31 
32 //
33 // constructors and destructor
34 //
36 {
37 
38  category_ = "Muon|RecoMuon|CosmicMuon|CosmicMuonLinksProducer";
39 
40  ParameterSet serviceParameters = iConfig.getParameter<ParameterSet>("ServiceParameters");
41 
42  theService = new MuonServiceProxy(serviceParameters);
43 
44  std::vector<edm::ParameterSet> theMapPSets = iConfig.getParameter<std::vector<edm::ParameterSet> >("Maps");
45  for (std::vector<edm::ParameterSet>::const_iterator iMPS = theMapPSets.begin();
46  iMPS != theMapPSets.end(); iMPS++) {
47 
48  edm::InputTag sTag = (*iMPS).getParameter<edm::InputTag>("subTrack");
49  edm::InputTag pTag = (*iMPS).getParameter<edm::InputTag>("parentTrack");
50 
51  edm::EDGetTokenT<reco::TrackCollection> subTrackTag = consumes<reco::TrackCollection>(sTag );
52  edm::EDGetTokenT<reco::TrackCollection> parentTrackTag = consumes<reco::TrackCollection>(pTag);
53  theTrackLinks.push_back( make_pair(subTrackTag, parentTrackTag) );
54  theTrackLinkNames.push_back( make_pair(sTag.label(), pTag.label()) );
55 
56 
57  LogDebug(category_) << "preparing map between " << sTag<<" & "<< pTag;
58  std::string mapname = sTag.label() + "To" + pTag.label();
59  produces<reco::TrackToTrackMap>(mapname);
60 
61 
62  }
63 
64 }
65 
67 {
68  if (theService) delete theService;
69 }
70 
71 
72 // ------------ method called to produce the data ------------
73 void
75 {
76  LogInfo(category_) << "Processing event number: " << iEvent.id();
77 
78  theService->update(iSetup);
79 
80  unsigned int counter= 0;
81  for(std::vector<std::pair<edm::EDGetTokenT<reco::TrackCollection>,edm::EDGetTokenT<reco::TrackCollection> > >::const_iterator iLink = theTrackLinks.begin();
82  iLink != theTrackLinks.end(); iLink++){
83  LogDebug(category_) << "making map between " << (*iLink).first<<" and "<< (*iLink).second;
84  std::string mapname = theTrackLinkNames[counter].first + "To" + theTrackLinkNames[counter].second;
86 
88  Handle<reco::TrackCollection> parentTracks;
89 
90  iEvent.getByToken((*iLink).first, subTracks);
91  iEvent.getByToken((*iLink).second, parentTracks);
92 
93  ttmap = mapTracks(subTracks, parentTracks);
94  LogTrace(category_) << "Mapped: "<<
95  theTrackLinkNames[counter].first <<" "<<subTracks->size()<< " and "<<theTrackLinkNames[counter].second<<" "<<parentTracks->size()<<", results: "<< ttmap.size() <<endl;
96 
97  auto_ptr<reco::TrackToTrackMap> trackToTrackmap(new reco::TrackToTrackMap(ttmap));
98  iEvent.put(trackToTrackmap, mapname);
99 
100  counter++;
101  }
102 
103 }
104 
107  for ( unsigned int position1 = 0; position1 != subTracks->size(); ++position1) {
108  TrackRef track1(subTracks, position1);
109  for ( unsigned int position2 = 0; position2 != parentTracks->size(); ++position2) {
110  TrackRef track2(parentTracks, position2);
111  int shared = sharedHits(*track1, *track2);
112  LogTrace(category_)<<"sharedHits "<<shared<<" track1 "<<track1->found()<<" track2 "<<track2->found()<<endl;
113 
114  if (shared > (track1->found())/2 ) map.insert(track1, track2);
115  }
116  }
117 
118  return map;
119 }
120 
121 int CosmicMuonLinksProducer::sharedHits(const reco::Track& track1, const reco::Track& track2) const {
122 
123  int match = 0;
124 
125  for (trackingRecHit_iterator hit1 = track1.recHitsBegin(); hit1 != track1.recHitsEnd(); ++hit1) {
126  if ( !(*hit1)->isValid() ) continue;
127  DetId id1 = (*hit1)->geographicalId();
128  if ( id1.det() != DetId::Muon ) continue; //ONLY MUON
129  LogTrace(category_)<<"first ID "<<id1.rawId()<<" "<<(*hit1)->localPosition()<<endl;
130  GlobalPoint pos1 = theService->trackingGeometry()->idToDet(id1)->surface().toGlobal((*hit1)->localPosition());
131 
132  for (trackingRecHit_iterator hit2 = track2.recHitsBegin(); hit2 != track2.recHitsEnd(); ++hit2) {
133 
134  if ( !(*hit2)->isValid() ) continue;
135 
136  DetId id2 = (*hit2)->geographicalId();
137  if ( id2.det() != DetId::Muon ) continue; //ONLY MUON
138 
139 // LogTrace(category_)<<"second ID "<<id2.rawId()<< (*hit2)->localPosition()<<endl;
140 
141  if (id2.rawId() != id1.rawId() ) continue;
142 
143  GlobalPoint pos2 = theService->trackingGeometry()->idToDet(id2)->surface().toGlobal((*hit2)->localPosition());
144  if ( ( pos1 - pos2 ).mag()< 10e-5 ) match++;
145 
146  }
147 
148  }
149 
150  return match;
151 
152 }
#define LogDebug(id)
T getParameter(std::string const &) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:446
virtual void produce(edm::Event &, const edm::EventSetup &)
int sharedHits(const reco::Track &track1, const reco::Track &track2) const
int sharedHits(const reco::GsfTrackRef &, const reco::GsfTrackRef &)
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:104
#define LogTrace(id)
Definition: DetId.h:18
size_type size() const
map size
void insert(const key_type &k, const data_type &v)
insert an association
std::string const & label() const
Definition: InputTag.h:42
edm::EventID id() const
Definition: EventBase.h:56
static std::atomic< unsigned int > counter
reco::TrackToTrackMap mapTracks(const edm::Handle< reco::TrackCollection > &, const edm::Handle< reco::TrackCollection > &) const
Detector det() const
get the detector field from this detid
Definition: DetId.h:35
CosmicMuonLinksProducer(const edm::ParameterSet &)
TrackingRecHitCollection::base::const_iterator trackingRecHit_iterator
iterator over a vector of reference to TrackingRecHit in the same collection
std::string match(BranchDescription const &a, BranchDescription const &b, std::string const &fileName)
trackingRecHit_iterator recHitsEnd() const
Iterator to last hit on the track.
Definition: Track.h:109