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