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,
17  bool
18 #ifdef EDM_ML_DEBUG
19  debug
20 #endif
21  ) {
22 
23  edm::SimTrackContainer::const_iterator itr = SimTk->end();
24  ;
25 
26  //Get the vector of PsimHits associated to TrackerRecHits and select the
27  //matching SimTrack on the basis of maximum occurance of trackIds
28  std::vector<unsigned int> trkId, trkOcc;
29  for (auto const& trkHit : pTrack->recHits()) {
30  std::vector<PSimHit> matchedSimIds = associate.associateHit(*trkHit);
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) {
42  trkId.push_back(tkId);
43  trkOcc.push_back(1);
44  }
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) {
64  maxTrkOcc = trkOcc[j];
65  idxMax = j;
66  }
67  }
68  matchSimTrk = trkId[idxMax];
69  for (auto simTrkItr = SimTk->begin(); simTrkItr != SimTk->end(); simTrkItr++) {
70  if (simTrkItr->trackId() == matchSimTrk) {
71 #ifdef EDM_ML_DEBUG
72  matchedId = simTrkItr->type();
73  if (debug)
74  std::cout << "matched trackId (maximum occurance) " << matchSimTrk << " type " << matchedId << std::endl;
75 #endif
76  itr = simTrkItr;
77  break;
78  }
79  }
80  }
81 
82 #ifdef EDM_ML_DEBUG
83  if (matchedId == 0 && debug) {
84  std::cout << "Could not find matched SimTrk and track history now " << std::endl;
85  }
86 #endif
87  return itr;
88  }
89 
90  //Returns a vector of TrackId originating from the matching SimTrack
91  std::vector<int> matchedSimTrackId(const edm::Event& iEvent,
94  const reco::Track* pTrack,
95  TrackerHitAssociator& associate,
96  bool debug) {
97  // get the matching SimTrack
98  edm::SimTrackContainer::const_iterator trkInfo =
99  spr::matchedSimTrack(iEvent, SimTk, SimVtx, pTrack, associate, debug);
100  unsigned int matchSimTrk = trkInfo->trackId();
101 #ifdef EDM_ML_DEBUG
102  if (debug)
103  std::cout << "matchedSimTrackId finds the SimTrk ID of the current track to be " << matchSimTrk << std::endl;
104 #endif
105  std::vector<int> matchTkid;
106  if (trkInfo->type() != 0) {
107  for (auto simTrkItr = SimTk->begin(); simTrkItr != SimTk->end(); simTrkItr++) {
108  if (validSimTrack(matchSimTrk, simTrkItr, SimTk, SimVtx, false))
109  matchTkid.push_back((int)simTrkItr->trackId());
110  }
111  }
112  return matchTkid;
113  }
114 
115  spr::simTkInfo matchedSimTrackInfo(unsigned int simTkId,
118  bool debug) {
120  for (edm::SimTrackContainer::const_iterator simTrkItr = SimTk->begin(); simTrkItr != SimTk->end(); simTrkItr++) {
121  if (simTkId == simTrkItr->trackId()) {
122  if (spr::validSimTrack(simTkId, simTrkItr, SimTk, SimVtx, debug)) {
123  info.found = true;
124  info.pdgId = simTrkItr->type();
125  info.charge = simTrkItr->charge();
126  } else {
127  edm::SimTrackContainer::const_iterator parentItr = spr::parentSimTrack(simTrkItr, SimTk, SimVtx, debug);
128 #ifdef EDM_ML_DEBUG
129  if (debug) {
130  if (parentItr != SimTk->end())
131  std::cout << "original parent of " << simTrkItr->trackId() << " " << parentItr->trackId() << ", "
132  << parentItr->type() << std::endl;
133  else
134  std::cout << "original parent of " << simTrkItr->trackId() << " not found" << std::endl;
135  }
136 #endif
137  if (parentItr != SimTk->end()) {
138  info.found = true;
139  info.pdgId = parentItr->type();
140  info.charge = parentItr->charge();
141  }
142  }
143  break;
144  }
145  }
146  return info;
147  }
148 
149  // Checks if this SimTrack=thisTrkItr originates from the one with trackId=simTkId
150  bool validSimTrack(unsigned int simTkId,
151  edm::SimTrackContainer::const_iterator thisTrkItr,
154  bool debug) {
155 #ifdef EDM_ML_DEBUG
156  if (debug)
157  std::cout << "Inside validSimTrack: trackId " << thisTrkItr->trackId() << " vtxIndex " << thisTrkItr->vertIndex()
158  << " to be matched to " << simTkId << std::endl;
159 #endif
160  //This track originates from simTkId
161  if (thisTrkItr->trackId() == simTkId)
162  return true;
163 
164  //Otherwise trace back the history using SimTracks and SimVertices
165  int vertIndex = thisTrkItr->vertIndex();
166  if (vertIndex == -1 || vertIndex >= (int)SimVtx->size())
167  return false;
168 
169  edm::SimVertexContainer::const_iterator simVtxItr = SimVtx->begin();
170  for (int iv = 0; iv < vertIndex; iv++)
171  simVtxItr++;
172  int parent = simVtxItr->parentIndex();
173 #ifdef EDM_ML_DEBUG
174  if (debug)
175  std::cout << "validSimTrack:: parent index " << parent << " ";
176 #endif
177  if (parent < 0 && simVtxItr != SimVtx->begin()) {
178  const math::XYZTLorentzVectorD pos1 = simVtxItr->position();
179  for (simVtxItr = SimVtx->begin(); simVtxItr != SimVtx->end(); ++simVtxItr) {
180  if (simVtxItr->parentIndex() > 0) {
181  const math::XYZTLorentzVectorD pos2 = pos1 - simVtxItr->position();
182  double dist = pos2.P();
183  if (dist < 0.001) {
184  parent = simVtxItr->parentIndex();
185  break;
186  }
187  }
188  }
189  }
190 #ifdef EDM_ML_DEBUG
191  if (debug)
192  std::cout << "final index " << parent << std::endl;
193  ;
194 #endif
195  for (edm::SimTrackContainer::const_iterator simTrkItr = SimTk->begin(); simTrkItr != SimTk->end(); simTrkItr++) {
196  if ((int)simTrkItr->trackId() == parent && simTrkItr != thisTrkItr)
197  return validSimTrack(simTkId, simTrkItr, SimTk, SimVtx, debug);
198  }
199 
200  return false;
201  }
202 
203  //Returns the parent of a SimTrack
204  edm::SimTrackContainer::const_iterator parentSimTrack(edm::SimTrackContainer::const_iterator thisTrkItr,
207  bool debug) {
208  edm::SimTrackContainer::const_iterator itr = SimTk->end();
209 
210  int vertIndex = thisTrkItr->vertIndex();
211 #ifdef EDM_ML_DEBUG
212  if (debug)
213  std::cout << "SimTrackParent " << thisTrkItr->trackId() << " Vertex " << vertIndex << " Type "
214  << thisTrkItr->type() << " Charge " << (int)thisTrkItr->charge() << std::endl;
215 #endif
216  if (vertIndex == -1)
217  return thisTrkItr;
218  else if (vertIndex >= (int)SimVtx->size())
219  return itr;
220 
221  edm::SimVertexContainer::const_iterator simVtxItr = SimVtx->begin();
222  for (int iv = 0; iv < vertIndex; iv++)
223  simVtxItr++;
224  int parent = simVtxItr->parentIndex();
225 
226  if (parent < 0 && simVtxItr != SimVtx->begin()) {
227  const math::XYZTLorentzVectorD pos1 = simVtxItr->position();
228  for (simVtxItr = SimVtx->begin(); simVtxItr != SimVtx->end(); ++simVtxItr) {
229  if (simVtxItr->parentIndex() > 0) {
230  const math::XYZTLorentzVectorD pos2 = pos1 - simVtxItr->position();
231  double dist = pos2.P();
232  if (dist < 0.001) {
233  parent = simVtxItr->parentIndex();
234  break;
235  }
236  }
237  }
238  }
239  for (edm::SimTrackContainer::const_iterator simTrkItr = SimTk->begin(); simTrkItr != SimTk->end(); simTrkItr++) {
240  if ((int)simTrkItr->trackId() == parent && simTrkItr != thisTrkItr)
241  return parentSimTrack(simTrkItr, SimTk, SimVtx, debug);
242  }
243 
244  return thisTrkItr;
245  }
246 
247 } // namespace spr
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:85
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)