CMS 3D CMS Logo

SimTrackManager.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: Application
4 // Class : SimTrackManager
5 //
6 // Implementation:
7 // <Notes on implementation>
8 //
9 // Original Author:
10 // Created: Fri Nov 25 17:44:19 EST 2005
11 //
12 
13 // system include files
14 #include <iostream>
15 
16 // user include files
21 
23 
24 #include "G4VProcess.hh"
25 
26 //#define DebugLog
27 
28 SimTrackManager::SimTrackManager(bool iCollapsePrimaryVertices)
29  : m_trksForThisEvent(nullptr),
30  m_nVertices(0),
31  m_collapsePrimaryVertices(iCollapsePrimaryVertices),
32  lastTrack(0),
33  lastHist(0),
34  theLHCTlink(nullptr) {}
35 
37  if (m_trksForThisEvent != nullptr)
38  deleteTracks();
39 }
40 
41 //
42 // member functions
43 //
45  if (m_trksForThisEvent == nullptr) {
47  } else {
48  for (unsigned int i = 0; i < m_trksForThisEvent->size(); i++) {
49  delete (*m_trksForThisEvent)[i];
50  }
51  delete m_trksForThisEvent;
53  }
56  std::vector<std::pair<int, int> >().swap(idsave);
57  ancestorList.clear();
58  lastTrack = 0;
59  lastHist = 0;
60 }
61 
63  for (unsigned int i = 0; i < m_trksForThisEvent->size(); i++) {
64  delete (*m_trksForThisEvent)[i];
65  }
66  delete m_trksForThisEvent;
67  m_trksForThisEvent = nullptr;
68 }
69 
72  using namespace std;
73  TrackWithHistory* trkH = trkWHist;
74  if (trkH == nullptr) {
75  edm::LogError("SimTrackManager") << " SimTrackManager::saveTrackAndItsBranch got 0 pointer ";
76  throw cms::Exception("SimTrackManager::saveTrackAndItsBranch") << " cannot handle hits for tracking";
77  }
78  trkH->save();
79  unsigned int parent = trkH->parentID();
80 
81  TrackContainer::const_iterator tk_itr = std::lower_bound(
82  (*m_trksForThisEvent).begin(), (*m_trksForThisEvent).end(), parent, SimTrackManager::StrictWeakOrdering());
83 
84  TrackWithHistory* tempTk = *tk_itr;
85 
86  if (tk_itr != m_trksForThisEvent->end() && (*tk_itr)->trackID() == parent) {
87  saveTrackAndItsBranch(tempTk);
88  }
89 }
90 
93 
94  // fill the map with the final mother-daughter relationship
95  idsave.swap(ancestorList);
96  stable_sort(idsave.begin(), idsave.end());
97 
98  std::vector<std::pair<int, int> >().swap(ancestorList);
99 
100  // to get a backward compatible order
101  stable_sort(m_trksForThisEvent->begin(), m_trksForThisEvent->end(), trkIDLess());
102 
103  // to reset the GenParticle ID of a SimTrack to its pre-LHCTransport value
104  resetGenID();
105 
106  reallyStoreTracks(simEvent);
107 }
108 
110  // loop over the (now ordered) vector and really save the tracks
111 #ifdef DebugLog
112  LogDebug("SimTrackManager") << "Inside the reallyStoreTracks method object to be stored = "
113  << m_trksForThisEvent->size();
114 #endif
115 
116  for (unsigned int it = 0; it < m_trksForThisEvent->size(); it++) {
117  TrackWithHistory* trkH = (*m_trksForThisEvent)[it];
118  // at this stage there is one vertex per track,
119  // so the vertex id of track N is also N
120  int ivertex = -1;
121  int ig;
122 
123  math::XYZVectorD pm(0., 0., 0.);
124  unsigned int iParentID = trkH->parentID();
125  for (unsigned int iit = 0; iit < m_trksForThisEvent->size(); ++iit) {
126  if ((*m_trksForThisEvent)[iit]->trackID() == iParentID) {
127  pm = (*m_trksForThisEvent)[iit]->momentum();
128  break;
129  }
130  }
131  ig = trkH->genParticleID();
132  ivertex = getOrCreateVertex(trkH, iParentID, simEvent);
133  std::map<uint32_t, std::pair<math::XYZVectorD, math::XYZTLorentzVectorD> >::const_iterator cit =
134  mapTkCaloStateInfo.find(trkH->trackID());
135  std::pair<math::XYZVectorD, math::XYZTLorentzVectorD> tcinfo;
136  if (cit != mapTkCaloStateInfo.end()) {
137  tcinfo = cit->second;
138  }
139  simEvent->add(new G4SimTrack(trkH->trackID(),
140  trkH->particleID(),
141  trkH->momentum(),
142  trkH->totalEnergy(),
143  ivertex,
144  ig,
145  pm,
146  tcinfo.first,
147  tcinfo.second));
148  }
149 }
150 
151 int SimTrackManager::getOrCreateVertex(TrackWithHistory* trkH, int iParentID, G4SimEvent* simEvent) {
152  int parent = iParentID;
153  int check = -1;
154 
155  for (std::vector<TrackWithHistory*>::const_iterator it = (*m_trksForThisEvent).begin();
156  it != (*m_trksForThisEvent).end();
157  it++) {
158  if ((*it)->trackID() == uint32_t(parent)) {
159  check = 0;
160  break;
161  }
162  }
163 
164  if (check == -1) {
165  parent = -1;
166  }
167 
168  VertexMap::const_iterator iterator = m_vertexMap.find(parent);
169  if (iterator != m_vertexMap.end()) {
170  // loop over saved vertices
171  for (unsigned int k = 0; k < m_vertexMap[parent].size(); k++) {
172  if (sqrt((trkH->vertexPosition() - (((m_vertexMap[parent])[k]).second)).Mag2()) < 0.001) {
173  return (((m_vertexMap[parent])[k]).first);
174  }
175  }
176  }
177 
178  unsigned int ptype = 0;
179  const G4VProcess* pr = trkH->creatorProcess();
180  if (pr) {
181  ptype = pr->GetProcessSubType();
182  }
183  simEvent->add(new G4SimVertex(trkH->vertexPosition(), trkH->globalTime(), parent, ptype));
185  m_nVertices++;
186  return (m_nVertices - 1);
187 }
188 
190  m_vertexMap.clear();
192  m_nVertices = 0;
193 }
194 
196  mapTkCaloStateInfo.clear();
197  std::map<uint32_t, std::pair<math::XYZVectorD, math::XYZTLorentzVectorD> >().swap(mapTkCaloStateInfo);
198 }
199 
200 int SimTrackManager::idSavedTrack(int id) const {
201  int idMother = id;
202  if (id > 0) {
203  unsigned int n = idsave.size();
204  if (0 < n) {
205  int jmax = n - 1;
206  int j, id1;
207 
208  // first loop forward
209  bool notFound = true;
210  for (j = 0; j <= jmax; ++j) {
211  if ((idsave[j]).first == idMother) {
212  id1 = (idsave[j]).second;
213  if (0 == id1 || id1 == idMother) {
214  return id1;
215  }
216  jmax = j - 1;
217  idMother = id1;
218  notFound = false;
219  break;
220  }
221  }
222  if (notFound) {
223  return 0;
224  }
225 
226  // recursive loop
227  do {
228  notFound = true;
229  // search ID scan backward
230  for (j = jmax; j >= 0; --j) {
231  if ((idsave[j]).first == idMother) {
232  id1 = (idsave[j]).second;
233  if (0 == id1 || id1 == idMother) {
234  return id1;
235  }
236  jmax = j - 1;
237  idMother = id1;
238  notFound = false;
239  break;
240  }
241  }
242  if (notFound) {
243  // ID not in the list of saved track - look into ancestors
244  jmax = ancestorList.size() - 1;
245  for (j = jmax; j >= 0; --j) {
246  if ((ancestorList[j]).first == idMother) {
247  idMother = (ancestorList[j]).second;
248  return idMother;
249  }
250  }
251  return 0;
252  }
253  } while (!notFound);
254  }
255  }
256  return idMother;
257 }
258 
260  if (!ancestorList.empty() && lastHist > ancestorList.size()) {
261  lastHist = ancestorList.size();
262  edm::LogError("SimTrackManager") << " SimTrackManager::fillMotherList track index corrupted";
263  }
264  /*
265  std::cout << "### SimTrackManager::fillMotherList: "
266  << idsave.size() << " saved; ancestor: " << lastHist
267  << " " << ancestorList.size() << std::endl;
268  for (unsigned int i = 0; i< idsave.size(); ++i) {
269  std::cout << " ISV: Track ID = " << (idsave[i]).first
270  << " Mother ID = " << (idsave[i]).second << std::endl;
271  }
272  */
273  for (unsigned int n = lastHist; n < ancestorList.size(); n++) {
274  int theMotherId = idSavedTrack((ancestorList[n]).first);
275  ancestorList[n].second = theMotherId;
276  /*
277  std::cout << " ANC: Track ID = " << (ancestorList[n]).first
278  << " Mother ID = " << (ancestorList[n]).second << std::endl;
279  */
280 #ifdef DebugLog
281  LogDebug("SimTrackManager") << "Track ID = " << (ancestorList[n]).first
282  << " Mother ID = " << (ancestorList[n]).second;
283 #endif
284  }
285 
286  lastHist = ancestorList.size();
287 
288  idsave.clear();
289 }
290 
292  if ((*m_trksForThisEvent).empty() && idsave.empty()) {
293  return;
294  }
295 
296 #ifdef DebugLog
297  LogDebug("SimTrackManager") << "SimTrackManager::cleanTracksWithHistory has " << idsave.size()
298  << " mother-daughter relationships stored with lastTrack = " << lastTrack;
299 #endif
300 
301  if (lastTrack > 0 && lastTrack >= (*m_trksForThisEvent).size()) {
302  lastTrack = 0;
303  edm::LogError("SimTrackManager") << " SimTrackManager::cleanTracksWithHistory track index corrupted";
304  }
305 
306  stable_sort(m_trksForThisEvent->begin() + lastTrack, m_trksForThisEvent->end(), trkIDLess());
307 
308  stable_sort(idsave.begin(), idsave.end());
309 
310 #ifdef DebugLog
311  LogDebug("SimTrackManager") << " SimTrackManager::cleanTracksWithHistory knows " << m_trksForThisEvent->size()
312  << " tracks with history before branching";
313  for (unsigned int it = 0; it < (*m_trksForThisEvent).size(); it++) {
314  LogDebug("SimTrackManager") << " 1 - Track in position " << it << " G4 track number "
315  << (*m_trksForThisEvent)[it]->trackID() << " mother "
316  << (*m_trksForThisEvent)[it]->parentID() << " status "
317  << (*m_trksForThisEvent)[it]->saved();
318  }
319 #endif
320 
321  for (unsigned int it = lastTrack; it < m_trksForThisEvent->size(); it++) {
322  TrackWithHistory* t = (*m_trksForThisEvent)[it];
323  if (t->saved()) {
325  }
326  }
327  unsigned int num = lastTrack;
328  for (unsigned int it = lastTrack; it < m_trksForThisEvent->size(); it++) {
329  TrackWithHistory* t = (*m_trksForThisEvent)[it];
330  int g4ID = t->trackID();
331  if (t->saved() == true) {
332  if (it > num)
333  (*m_trksForThisEvent)[num] = t;
334  num++;
335  for (unsigned int itr = 0; itr < idsave.size(); itr++) {
336  if ((idsave[itr]).first == g4ID) {
337  (idsave[itr]).second = g4ID;
338  break;
339  }
340  }
341  } else {
342  delete t;
343  }
344  }
345 
346  (*m_trksForThisEvent).resize(num);
347 
348 #ifdef DebugLog
349  LogDebug("SimTrackManager") << " AFTER CLEANING, I GET " << (*m_trksForThisEvent).size()
350  << " tracks to be saved persistently";
351  for (unsigned int it = 0; it < (*m_trksForThisEvent).size(); it++) {
352  LogDebug("SimTrackManager") << " Track in position " << it << " G4 track number "
353  << (*m_trksForThisEvent)[it]->trackID() << " mother "
354  << (*m_trksForThisEvent)[it]->parentID() << " Status "
355  << (*m_trksForThisEvent)[it]->saved() << " id "
356  << (*m_trksForThisEvent)[it]->particleID()
357  << " E(MeV)= " << (*m_trksForThisEvent)[it]->totalEnergy();
358  }
359 #endif
360 
361  fillMotherList();
362 
363  lastTrack = (*m_trksForThisEvent).size();
364 }
365 
367  if (theLHCTlink == nullptr)
368  return;
369 
370  for (unsigned int it = 0; it < m_trksForThisEvent->size(); it++) {
371  TrackWithHistory* trkH = (*m_trksForThisEvent)[it];
372  int genParticleID_ = trkH->genParticleID();
373  if (genParticleID_ == -1) {
374  continue;
375  } else {
376  for (unsigned int itrlink = 0; itrlink < (*theLHCTlink).size(); itrlink++) {
377  if ((*theLHCTlink)[itrlink].afterHector() == genParticleID_) {
378  trkH->setGenParticleID((*theLHCTlink)[itrlink].beforeHector());
379  continue;
380  }
381  }
382  }
383  }
384 
385  theLHCTlink = nullptr;
386 }
#define LogDebug(id)
int getOrCreateVertex(TrackWithHistory *, int, G4SimEvent *simEvent)
bool saved() const
std::map< int, MapVertexPositionVector > MotherParticleToVertexMap
#define nullptr
double totalEnergy() const
SimTrackManager(bool iCollapsePrimaryVertices=false)
void setGenParticleID(int i)
void cleanTkCaloStateInfoMap()
void swap(Association< C > &lhs, Association< C > &rhs)
Definition: Association.h:116
unsigned int lastHist
virtual ~SimTrackManager()
TrackContainer * m_trksForThisEvent
U second(std::pair< T, U > const &p)
void add(G4SimTrack *t)
Definition: G4SimEvent.h:33
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double > > XYZVectorD
spatial vector with cartesian internal representation
Definition: Vector3D.h:8
std::vector< TrackWithHistory * > TrackContainer
Definition: TrackContainer.h:8
T sqrt(T t)
Definition: SSEVec.h:18
unsigned int trackID() const
void cleanTracksWithHistory()
MotherParticleToVertexMap m_vertexMap
int particleID() const
void storeTracks(G4SimEvent *simEvent)
int idSavedTrack(int) const
int k[5][pyjets_maxn]
int parentID() const
const G4VProcess * creatorProcess() const
std::vector< std::pair< int, int > > idsave
std::map< uint32_t, std::pair< math::XYZVectorD, math::XYZTLorentzVectorD > > mapTkCaloStateInfo
std::pair< int, math::XYZVectorD > MapVertexPosition
this map contains association between vertex number and position
const math::XYZVectorD & momentum() const
void saveTrackAndItsBranch(TrackWithHistory *)
this saves a track and all its parents looping over the non ordered vector
unsigned int lastTrack
double globalTime() const
static const GlobalPoint notFound(0, 0, 0)
const edm::LHCTransportLinkContainer * theLHCTlink
void reallyStoreTracks(G4SimEvent *simEvent)
std::vector< std::pair< int, int > > ancestorList
const math::XYZVectorD & vertexPosition() const
def check(config)
Definition: trackerTree.py:14
int genParticleID() const