CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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
20 
22 
23 #include "G4VProcess.hh"
24 
25 //#define DebugLog
26 
27 //
28 // constants, enums and typedefs
29 //
30 
31 //
32 // static data member definitions
33 //
34 
35 //
36 // constructors and destructor
37 //
38 SimTrackManager::SimTrackManager(bool iCollapsePrimaryVertices) :
39  m_trksForThisEvent(0),m_nVertices(0),
40  m_collapsePrimaryVertices(iCollapsePrimaryVertices),
41  lastTrack(0),lastHist(0),theLHCTlink(0){}
42 
43 
45 {
46  if ( m_trksForThisEvent != 0 ) deleteTracks() ;
47 }
48 
49 //
50 // member functions
51 //
53 {
55  else
56  {
57  for (unsigned int i = 0; i < m_trksForThisEvent->size(); i++) {
58  delete (*m_trksForThisEvent)[i];
59  }
60  delete m_trksForThisEvent;
62  }
65  std::vector<std::pair <int, int> >().swap(idsave);
66  ancestorList.clear();
67  lastTrack=0;
68  lastHist=0;
69 }
70 
72 {
73  for (unsigned int i = 0; i < m_trksForThisEvent->size(); i++) {
74  delete (*m_trksForThisEvent)[i];
75  }
76  delete m_trksForThisEvent;
78 }
79 
82 {
83  using namespace std;
84  TrackWithHistory * trkH = trkWHist;
85  if (trkH == 0)
86  {
87  edm::LogError("SimTrackManager")
88  << " SimTrackManager::saveTrackAndItsBranch got 0 pointer ";
89  abort();
90  }
91  trkH->save();
92  unsigned int parent = trkH->parentID();
93 
94  TrackContainer::const_iterator tk_itr =
95  std::lower_bound((*m_trksForThisEvent).begin(),(*m_trksForThisEvent).end(),
97 
98  TrackWithHistory * tempTk = *tk_itr;
99 
100  if (tk_itr!=m_trksForThisEvent->end() && (*tk_itr)->trackID()==parent) {
101  saveTrackAndItsBranch(tempTk);
102  }
103 }
104 
106 {
108 
109  // fill the map with the final mother-daughter relationship
110  idsave.swap(ancestorList);
111  stable_sort(idsave.begin(),idsave.end());
112 
113  std::vector<std::pair<int,int> >().swap(ancestorList);
114 
115  // to get a backward compatible order
116  stable_sort(m_trksForThisEvent->begin(),m_trksForThisEvent->end(),trkIDLess());
117 
118  // to reset the GenParticle ID of a SimTrack to its pre-LHCTransport value
119  resetGenID();
120 
121  reallyStoreTracks(simEvent);
122 }
123 
125 {
126  // loop over the (now ordered) vector and really save the tracks
127 #ifdef DebugLog
128  LogDebug("SimTrackManager")
129  << "Inside the reallyStoreTracks method object to be stored = "
130  << m_trksForThisEvent->size();
131 #endif
132 
133  for (unsigned int it = 0; it < m_trksForThisEvent->size(); it++)
134  {
135  TrackWithHistory * trkH = (*m_trksForThisEvent)[it];
136  // at this stage there is one vertex per track,
137  // so the vertex id of track N is also N
138  int ivertex = -1;
139  int ig;
140 
141  math::XYZVectorD pm(0.,0.,0.);
142  unsigned int iParentID = trkH->parentID();
143  for(unsigned int iit = 0; iit < m_trksForThisEvent->size(); iit++)
144  {
145  if((*m_trksForThisEvent)[iit]->trackID()==iParentID){
146  pm = (*m_trksForThisEvent)[iit]->momentum();
147  break;
148  }
149  }
150  ig = trkH->genParticleID();
151  ivertex = getOrCreateVertex(trkH,iParentID,simEvent);
152  std::map<uint32_t,std::pair<math::XYZVectorD,math::XYZTLorentzVectorD> >::const_iterator cit =
153  mapTkCaloStateInfo.find(trkH->trackID());
154  std::pair<math::XYZVectorD,math::XYZTLorentzVectorD> tcinfo;
155  if (cit != mapTkCaloStateInfo.end()){
156  tcinfo = cit->second;
157  }
158  simEvent->add(new G4SimTrack(trkH->trackID(),trkH->particleID(),
159  trkH->momentum(),trkH->totalEnergy(),
160  ivertex,ig,pm,tcinfo.first,tcinfo.second));
161  }
162 }
163 
165  G4SimEvent * simEvent)
166 {
167  int parent = iParentID;
168  int check = -1;
169 
170  for( std::vector<TrackWithHistory*>::const_iterator it = (*m_trksForThisEvent).begin();
171  it!= (*m_trksForThisEvent).end();it++)
172  {
173  if ((*it)->trackID() == uint32_t(parent)) {
174  check = 0;
175  break;
176  }
177  }
178 
179  if(check==-1) { parent = -1; }
180 
181  VertexMap::const_iterator iterator = m_vertexMap.find(parent);
182  if (iterator != m_vertexMap.end()) {
183 
184  // loop over saved vertices
185  for(unsigned int k=0; k<m_vertexMap[parent].size(); k++){
186  if(sqrt((trkH->vertexPosition()-(((m_vertexMap[parent])[k]).second)).Mag2())<0.001)
187  { return (((m_vertexMap[parent])[k]).first); }
188  }
189  }
190 
191  unsigned int ptype = 0;
192  const G4VProcess* pr = trkH->creatorProcess();
193  if(pr) { ptype = pr->GetProcessSubType(); }
194  simEvent->add(new G4SimVertex(trkH->vertexPosition(),trkH->globalTime(),parent,ptype));
196  m_nVertices++;
197  return (m_nVertices-1);
198 
199 }
200 
202 {
203  m_vertexMap.clear();
205  m_nVertices=0;
206 }
207 
209 {
210  mapTkCaloStateInfo.clear();
211  std::map<uint32_t,std::pair<math::XYZVectorD,math::XYZTLorentzVectorD > >().swap(mapTkCaloStateInfo);
212 }
213 
215 {
216  int idMother = id;
217  if(id > 0) {
218  unsigned int n = idsave.size();
219  if(0 < n) {
220  int jmax = n - 1;
221  int j, id1;
222 
223  // first loop forward
224  bool notFound = true;
225  for(j=0; j<=jmax; ++j) {
226  if((idsave[j]).first == idMother) {
227  id1 = (idsave[j]).second;
228  if(0 == id1 || id1 == idMother) { return id1; }
229  jmax = j - 1;
230  idMother = id1;
231  notFound = false;
232  break;
233  }
234  }
235  if(notFound) { return 0; }
236 
237  // recursive loop
238  do {
239 
240  notFound = true;
241  // search ID scan backward
242  for(j=jmax; j>=0; --j) {
243  if((idsave[j]).first == idMother) {
244  id1 = (idsave[j]).second;
245  if(0 == id1 || id1 == idMother) { return id1; }
246  jmax = j - 1;
247  idMother = id1;
248  notFound = false;
249  break;
250  }
251  }
252  if(notFound) {
253  // ID not in the list of saved track - look into ancestors
254  jmax = ancestorList.size()-1;
255  for(j=jmax; j>=0; --j) {
256  if((ancestorList[j]).first == idMother) {
257  idMother = (ancestorList[j]).second;
258  return idMother;
259  }
260  }
261  return 0;
262  }
263  } while (!notFound);
264  }
265  }
266  return idMother;
267 }
268 
270 {
271  if ( ancestorList.size() > 0 && lastHist > ancestorList.size() ) {
272  lastHist = ancestorList.size();
273  edm::LogError("SimTrackManager")
274  << " SimTrackManager::fillMotherList track index corrupted";
275  }
276  /*
277  std::cout << "### SimTrackManager::fillMotherList: "
278  << idsave.size() << " saved; ancestor: " << lastHist
279  << " " << ancestorList.size() << std::endl;
280  for (unsigned int i = 0; i< idsave.size(); ++i) {
281  std::cout << " ISV: Track ID = " << (idsave[i]).first
282  << " Mother ID = " << (idsave[i]).second << std::endl;
283  }
284  */
285  for (unsigned int n = lastHist; n < ancestorList.size(); n++) {
286 
287  int theMotherId = idSavedTrack((ancestorList[n]).first);
288  ancestorList[n].second = theMotherId;
289  /*
290  std::cout << " ANC: Track ID = " << (ancestorList[n]).first
291  << " Mother ID = " << (ancestorList[n]).second << std::endl;
292  */
293 #ifdef DebugLog
294  LogDebug("SimTrackManager") << "Track ID = " << (ancestorList[n]).first
295  << " Mother ID = " << (ancestorList[n]).second;
296 #endif
297  }
298 
299  lastHist = ancestorList.size();
300 
301  idsave.clear();
302 
303 }
304 
306 
307  using namespace std;
308 
309  if ((*m_trksForThisEvent).size() == 0 && idsave.size() == 0) { return; }
310 
311 #ifdef DebugLog
312  LogDebug("SimTrackManager")
313  << "SimTrackManager::cleanTracksWithHistory has "
314  << idsave.size()
315  << " mother-daughter relationships stored with lastTrack = " << lastTrack;
316 #endif
317 
318  if ( lastTrack > 0 && lastTrack >= (*m_trksForThisEvent).size() ) {
319  lastTrack = 0;
320  edm::LogError("SimTrackManager")
321  << " SimTrackManager::cleanTracksWithHistory track index corrupted";
322  }
323 
324  stable_sort(m_trksForThisEvent->begin()+lastTrack,m_trksForThisEvent->end(),trkIDLess());
325 
326  stable_sort(idsave.begin(),idsave.end());
327 
328 #ifdef DebugLog
329  LogDebug("SimTrackManager")
330  << " SimTrackManager::cleanTracksWithHistory knows " << m_trksForThisEvent->size()
331  << " tracks with history before branching";
332  for (unsigned int it =0; it <(*m_trksForThisEvent).size(); it++) {
333  LogDebug("SimTrackManager")
334  << " 1 - Track in position " << it << " G4 track number "
335  << (*m_trksForThisEvent)[it]->trackID()
336  << " mother " << (*m_trksForThisEvent)[it]->parentID()
337  << " status " << (*m_trksForThisEvent)[it]->saved();
338  }
339 #endif
340 
341  for (unsigned int it = lastTrack; it < m_trksForThisEvent->size(); it++)
342  {
343  TrackWithHistory * t = (*m_trksForThisEvent)[it];
344  if (t->saved()) { saveTrackAndItsBranch(t); }
345  }
346  unsigned int num = lastTrack;
347  for (unsigned int it = lastTrack; it < m_trksForThisEvent->size(); it++)
348  {
349  TrackWithHistory * t = (*m_trksForThisEvent)[it];
350  int g4ID = t->trackID();
351  if (t->saved() == true)
352  {
353  if (it>num) (*m_trksForThisEvent)[num] = t;
354  num++;
355  for (unsigned int itr=0; itr<idsave.size(); itr++) {
356  if ((idsave[itr]).first == g4ID) {
357  (idsave[itr]).second = g4ID;
358  break;
359  }
360  }
361  }
362  else
363  {
364  delete t;
365  }
366  }
367 
368  (*m_trksForThisEvent).resize(num);
369 
370 #ifdef DebugLog
371  LogDebug("SimTrackManager")
372  << " AFTER CLEANING, I GET " << (*m_trksForThisEvent).size()
373  << " tracks to be saved persistently";
374  for (unsigned int it = 0; it < (*m_trksForThisEvent).size(); it++) {
375  LogDebug("SimTrackManager")
376  << " Track in position " << it
377  << " G4 track number " << (*m_trksForThisEvent)[it]->trackID()
378  << " mother " << (*m_trksForThisEvent)[it]->parentID()
379  << " Status " << (*m_trksForThisEvent)[it]->saved();
380  }
381 #endif
382 
383  fillMotherList();
384 
385  lastTrack = (*m_trksForThisEvent).size();
386 
387 }
388 
390 {
391  if ( theLHCTlink == 0 ) return;
392 
393  for (unsigned int it = 0; it < m_trksForThisEvent->size(); it++)
394  {
395  TrackWithHistory * trkH = (*m_trksForThisEvent)[it];
396  int genParticleID_ = trkH->genParticleID();
397  if ( genParticleID_ == -1 ) { continue; }
398  else {
399  for ( unsigned int itrlink = 0; itrlink < (*theLHCTlink).size(); itrlink++ ) {
400  if ( (*theLHCTlink)[itrlink].afterHector() == genParticleID_ ) {
401  trkH->setGenParticleID( (*theLHCTlink)[itrlink].beforeHector() );
402  continue;
403  }
404  }
405  }
406  }
407 
408  theLHCTlink = 0;
409 
410 }
#define LogDebug(id)
void swap(ora::Record &rh, ora::Record &lh)
Definition: Record.h:70
int i
Definition: DBlmapReader.cc:9
int getOrCreateVertex(TrackWithHistory *, int, G4SimEvent *simEvent)
list parent
Definition: dbtoconf.py:74
bool saved() const
double totalEnergy() const
SimTrackManager(bool iCollapsePrimaryVertices=false)
void setGenParticleID(int i)
void cleanTkCaloStateInfoMap()
std::map< uint32_t, std::pair< math::XYZVectorD, math::XYZTLorentzVectorD > > mapTkCaloStateInfo
unsigned int lastHist
virtual ~SimTrackManager()
TrackContainer * m_trksForThisEvent
U second(std::pair< T, U > const &p)
void add(G4SimTrack *t)
Definition: G4SimEvent.h:35
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:48
bool check(const DataFrame &df, bool capcheck, bool dvercheck)
std::vector< std::pair< int, int > > idsave
unsigned int trackID() const
void cleanTracksWithHistory()
int j
Definition: DBlmapReader.cc:9
std::pair< int, math::XYZVectorD > MapVertexPosition
this map contains association between vertex number and position
MotherParticleToVertexMap m_vertexMap
int particleID() const
void storeTracks(G4SimEvent *simEvent)
int idSavedTrack(int) const
bool first
Definition: L1TdeRCT.cc:79
int k[5][pyjets_maxn]
int parentID() const
const G4VProcess * creatorProcess() const
std::map< int, MapVertexPositionVector > MotherParticleToVertexMap
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
int genParticleID() const