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