CMS 3D CMS Logo

MatchingSimTrack.cc
Go to the documentation of this file.
3 
5 
6 #include <iostream>
7 
8 //#define EDM_ML_DEBUG
9 
10 namespace spr{
11 
12  edm::SimTrackContainer::const_iterator matchedSimTrack(const edm::Event& iEvent,
15  const reco::Track* pTrack,
16  TrackerHitAssociator& associate, bool
17 #ifdef EDM_ML_DEBUG
18  debug
19 #endif
20  ) {
21 
22  edm::SimTrackContainer::const_iterator itr = SimTk->end();;
23 
24  edm::SimTrackContainer::const_iterator simTrkItr;
25  edm::SimVertexContainer::const_iterator simVtxItr;
26 
27  //Get the vector of PsimHits associated to TrackerRecHits and select the
28  //matching SimTrack on the basis of maximum occurance of trackIds
29  std::vector<unsigned int> trkId, trkOcc;
30  int i=0;
31  for (trackingRecHit_iterator iTrkHit = pTrack->recHitsBegin(); iTrkHit != pTrack->recHitsEnd(); ++iTrkHit, ++i) {
32 
33  std::vector<PSimHit> matchedSimIds = associate.associateHit((**iTrkHit));
34  for (unsigned int isim=0; isim<matchedSimIds.size(); isim++) {
35  unsigned tkId = matchedSimIds[isim].trackId();
36  bool found = false;
37  for (unsigned int j=0; j<trkId.size(); j++) {
38  if ( tkId == trkId[j] ) {
39  trkOcc[j]++;
40  found = true;
41  break;
42  }
43  }
44  if (!found) { trkId.push_back(tkId); trkOcc.push_back(1); }
45  }
46  }
47 
48 #ifdef EDM_ML_DEBUG
49  if (debug) {
50  std::cout << "Reconstructed Track with " << i << " recHits.";
51  for (unsigned int isim=0; isim<trkId.size(); isim++){
52  std::cout << "\n trkId " << trkId[isim] << " Occurance " << trkOcc[isim] << ", ";
53  }
54  std::cout << std::endl;
55  }
56  int matchedId=0;
57 #endif
58 
59  unsigned int matchSimTrk=0;
60  if (!trkOcc.empty()) {
61  unsigned int maxTrkOcc=0, idxMax=0;
62  for(unsigned int j=0; j<trkOcc.size(); j++) {
63  if(trkOcc[j] > maxTrkOcc ) { maxTrkOcc = trkOcc[j]; idxMax = j; }
64  }
65  matchSimTrk = trkId[idxMax];
66  for (simTrkItr = SimTk->begin(); simTrkItr!= SimTk->end(); simTrkItr++) {
67  if ( simTrkItr->trackId() == matchSimTrk ) {
68 #ifdef EDM_ML_DEBUG
69  matchedId = simTrkItr->type();
70  if (debug) std::cout << "matched trackId (maximum occurance) " << matchSimTrk << " type " << matchedId << std::endl;
71 #endif
72  itr = simTrkItr;
73  break;
74  }
75  }
76  }
77 
78 #ifdef EDM_ML_DEBUG
79  if (matchedId==0 && debug) {
80  std::cout << "Could not find matched SimTrk and track history now " << std::endl;
81  }
82 #endif
83  return itr;
84  }
85 
86  //Returns a vector of TrackId originating from the matching SimTrack
88 
89  // get the matching SimTrack
90  edm::SimTrackContainer::const_iterator trkInfo = spr::matchedSimTrack(iEvent, SimTk, SimVtx, pTrack,associate, debug);
91  unsigned int matchSimTrk = trkInfo->trackId();
92 #ifdef EDM_ML_DEBUG
93  if (debug) std::cout << "matchedSimTrackId finds the SimTrk ID of the current track to be " << matchSimTrk << std::endl;
94 #endif
95  std::vector<int> matchTkid;
96  if( trkInfo->type() != 0) {
97  edm::SimTrackContainer::const_iterator simTrkItr;
98  for(simTrkItr = SimTk->begin(); simTrkItr!= SimTk->end(); simTrkItr++){
99  if (validSimTrack(matchSimTrk, simTrkItr, SimTk, SimVtx, false))
100  matchTkid.push_back((int)simTrkItr->trackId());
101  }
102  }
103  return matchTkid;
104  }
105 
108  for (edm::SimTrackContainer::const_iterator simTrkItr=SimTk->begin();
109  simTrkItr!= SimTk->end(); simTrkItr++) {
110  if (simTkId == simTrkItr->trackId()) {
111  if (spr::validSimTrack(simTkId, simTrkItr, SimTk, SimVtx, debug)) {
112  info.found = true;
113  info.pdgId = simTrkItr->type();
114  info.charge = simTrkItr->charge();
115  } else {
116  edm::SimTrackContainer::const_iterator parentItr = spr::parentSimTrack(simTrkItr, SimTk, SimVtx, debug);
117 #ifdef EDM_ML_DEBUG
118  if (debug) {
119  if (parentItr != SimTk->end() ) std::cout << "original parent of " << simTrkItr->trackId() << " " << parentItr->trackId() << ", " << parentItr->type() << std::endl;
120  else std::cout << "original parent of " << simTrkItr->trackId() << " not found" << std::endl;
121  }
122 #endif
123  if (parentItr != SimTk->end()) {
124  info.found = true;
125  info.pdgId = parentItr->type();
126  info.charge = parentItr->charge();
127  }
128  }
129  break;
130  }
131  }
132  return info;
133  }
134 
135  // Checks if this SimTrack=thisTrkItr originates from the one with trackId=simTkId
136  bool validSimTrack(unsigned int simTkId, edm::SimTrackContainer::const_iterator thisTrkItr, edm::Handle<edm::SimTrackContainer>& SimTk, edm::Handle<edm::SimVertexContainer>& SimVtx, bool debug){
137 
138 #ifdef EDM_ML_DEBUG
139  if (debug) std::cout << "Inside validSimTrack: trackId " << thisTrkItr->trackId() << " vtxIndex " << thisTrkItr->vertIndex() << " to be matched to " << simTkId << std::endl;
140 #endif
141  //This track originates from simTkId
142  if (thisTrkItr->trackId() == simTkId) return true;
143 
144  //Otherwise trace back the history using SimTracks and SimVertices
145  int vertIndex = thisTrkItr->vertIndex();
146  if (vertIndex == -1 || vertIndex >= (int)SimVtx->size()) return false;
147 
148  edm::SimVertexContainer::const_iterator simVtxItr= SimVtx->begin();
149  for (int iv=0; iv<vertIndex; iv++) simVtxItr++;
150  int parent = simVtxItr->parentIndex();
151 #ifdef EDM_ML_DEBUG
152  if (debug) std::cout << "validSimTrack:: parent index " << parent <<" ";
153 #endif
154  if (parent < 0 && simVtxItr != SimVtx->begin()) {
155  const math::XYZTLorentzVectorD pos1 = simVtxItr->position();
156  for (simVtxItr=SimVtx->begin(); simVtxItr!=SimVtx->end(); ++simVtxItr) {
157  if (simVtxItr->parentIndex() > 0) {
158  const math::XYZTLorentzVectorD pos2 = pos1 - simVtxItr->position();
159  double dist = pos2.P();
160  if (dist < 0.001) {
161  parent = simVtxItr->parentIndex();
162  break;
163  }
164  }
165  }
166  }
167 #ifdef EDM_ML_DEBUG
168  if (debug) std::cout << "final index " << parent << std::endl;;
169 #endif
170  for(edm::SimTrackContainer::const_iterator simTrkItr= SimTk->begin(); simTrkItr!= SimTk->end(); simTrkItr++){
171  if ((int)simTrkItr->trackId() == parent && simTrkItr != thisTrkItr) return validSimTrack(simTkId, simTrkItr, SimTk, SimVtx, debug) ;
172  }
173 
174  return false;
175 
176  }
177 
178  //Returns the parent of a SimTrack
179  edm::SimTrackContainer::const_iterator parentSimTrack(edm::SimTrackContainer::const_iterator thisTrkItr, edm::Handle<edm::SimTrackContainer>& SimTk, edm::Handle<edm::SimVertexContainer>& SimVtx, bool debug){
180 
181  edm::SimTrackContainer::const_iterator itr = SimTk->end();
182 
183  int vertIndex = thisTrkItr->vertIndex();
184 #ifdef EDM_ML_DEBUG
185  if (debug) std::cout << "SimTrackParent " << thisTrkItr->trackId() << " Vertex " << vertIndex << " Type " << thisTrkItr->type() << " Charge " << (int)thisTrkItr->charge() << std::endl;
186 #endif
187  if( vertIndex == -1 ) return thisTrkItr;
188  else if (vertIndex >= (int)SimVtx->size()) return itr;
189 
190  edm::SimVertexContainer::const_iterator simVtxItr= SimVtx->begin();
191  for (int iv=0; iv<vertIndex; iv++) simVtxItr++;
192  int parent = simVtxItr->parentIndex();
193 
194  if (parent < 0 && simVtxItr != SimVtx->begin()) {
195  const math::XYZTLorentzVectorD pos1 = simVtxItr->position();
196  for (simVtxItr=SimVtx->begin(); simVtxItr!=SimVtx->end(); ++simVtxItr) {
197  if (simVtxItr->parentIndex() > 0) {
198  const math::XYZTLorentzVectorD pos2 = pos1 - simVtxItr->position();
199  double dist = pos2.P();
200  if (dist < 0.001) {
201  parent = simVtxItr->parentIndex();
202  break;
203  }
204  }
205  }
206  }
207  for (edm::SimTrackContainer::const_iterator simTrkItr= SimTk->begin(); simTrkItr!= SimTk->end(); simTrkItr++){
208  if ((int)simTrkItr->trackId() == parent && simTrkItr != thisTrkItr) return parentSimTrack(simTrkItr, SimTk, SimVtx, debug);
209  }
210 
211  return thisTrkItr;
212  }
213 
214 
215 }
static const TGPicture * info(bool iBackgroundIsBlack)
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:14
int iEvent
Definition: GenABIO.cc:230
edm::SimTrackContainer::const_iterator matchedSimTrack(const edm::Event &iEvent, edm::Handle< edm::SimTrackContainer > &SimTk, edm::Handle< edm::SimVertexContainer > &SimVtx, const reco::Track *pTrack, TrackerHitAssociator &associate, bool debug=false)
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:104
#define debug
Definition: HDRShower.cc:19
simTkInfo matchedSimTrackInfo(unsigned int simTkId, edm::Handle< edm::SimTrackContainer > &SimTk, edm::Handle< edm::SimVertexContainer > &SimVtx, bool debug=false)
#define begin
Definition: vmac.h:30
#define EDM_ML_DEBUG
std::vector< int > matchedSimTrackId(const edm::Event &, edm::Handle< edm::SimTrackContainer > &SimTk, edm::Handle< edm::SimVertexContainer > &SimVtx, const reco::Track *pTrack, TrackerHitAssociator &associate, bool debug=false)
std::vector< PSimHit > associateHit(const TrackingRecHit &thit) const
bool validSimTrack(unsigned int simTkId, edm::SimTrackContainer::const_iterator thisTrkItr, edm::Handle< edm::SimTrackContainer > &SimTk, edm::Handle< edm::SimVertexContainer > &SimVtx, bool debug=false)
TrackingRecHitCollection::base::const_iterator trackingRecHit_iterator
iterator over a vector of reference to TrackingRecHit in the same collection
edm::SimTrackContainer::const_iterator parentSimTrack(edm::SimTrackContainer::const_iterator thisTrkItr, edm::Handle< edm::SimTrackContainer > &SimTk, edm::Handle< edm::SimVertexContainer > &SimVtx, bool debug=false)
trackingRecHit_iterator recHitsEnd() const
Iterator to last hit on the track.
Definition: Track.h:109