CMS 3D CMS Logo

List of all members | Classes | Public Types | Public Member Functions | Private Member Functions | Private Attributes
SimTrackManager Class Reference

#include <SimG4Core/Notification/interface/SimTrackManager.h>

Classes

class  StrictWeakOrdering
 

Public Types

typedef std::pair< int, math::XYZVectorDMapVertexPosition
 this map contains association between vertex number and position More...
 
typedef std::vector< std::pair< int, math::XYZVectorD > > MapVertexPositionVector
 
typedef std::map< int, MapVertexPositionVectorMotherParticleToVertexMap
 
typedef MotherParticleToVertexMap VertexMap
 

Public Member Functions

void addTkCaloStateInfo (uint32_t t, const std::pair< math::XYZVectorD, math::XYZTLorentzVectorD > &p)
 
void addTrack (TrackWithHistory *iTrack, bool inHistory, bool withAncestor)
 
void cleanTkCaloStateInfoMap ()
 
void cleanTracksWithHistory ()
 
void deleteTracks ()
 
TrackWithHistorygetTrackByID (unsigned int trackID, bool strict=false) const
 
int giveMotherNeeded (int i) const
 
const SimTrackManageroperator= (const SimTrackManager &)=delete
 
void reset ()
 
void setCollapsePrimaryVertices (bool iSet)
 
void setLHCTransportLink (const edm::LHCTransportLinkContainer *thisLHCTlink)
 
 SimTrackManager (bool iCollapsePrimaryVertices=false)
 
 SimTrackManager (const SimTrackManager &)=delete
 
void storeTracks (G4SimEvent *simEvent)
 
const TrackContainertrackContainer () const
 
bool trackExists (unsigned int i) const
 
virtual ~SimTrackManager ()
 

Private Member Functions

void cleanVertexMap ()
 
void fillMotherList ()
 
int getOrCreateVertex (TrackWithHistory *, int, G4SimEvent *simEvent)
 
int idSavedTrack (int) const
 
void reallyStoreTracks (G4SimEvent *simEvent)
 
void resetGenID ()
 
void saveTrackAndItsBranch (TrackWithHistory *)
 this saves a track and all its parents looping over the non ordered vector More...
 

Private Attributes

std::vector< std::pair< int, int > > ancestorList
 
std::vector< std::pair< int, int > > idsave
 
unsigned int lastHist
 
unsigned int lastTrack
 
bool m_collapsePrimaryVertices
 
int m_nVertices
 
bool m_SaveSimTracks
 
TrackContainerm_trksForThisEvent
 
MotherParticleToVertexMap m_vertexMap
 
std::map< uint32_t, std::pair< math::XYZVectorD, math::XYZTLorentzVectorD > > mapTkCaloStateInfo
 
const edm::LHCTransportLinkContainertheLHCTlink
 

Detailed Description

Description: Holds tracking information used by the sensitive detectors

Usage: <usage>

Definition at line 35 of file SimTrackManager.h.

Member Typedef Documentation

◆ MapVertexPosition

this map contains association between vertex number and position

Definition at line 43 of file SimTrackManager.h.

◆ MapVertexPositionVector

typedef std::vector<std::pair<int, math::XYZVectorD> > SimTrackManager::MapVertexPositionVector

Definition at line 44 of file SimTrackManager.h.

◆ MotherParticleToVertexMap

Definition at line 45 of file SimTrackManager.h.

◆ VertexMap

Definition at line 46 of file SimTrackManager.h.

Constructor & Destructor Documentation

◆ SimTrackManager() [1/2]

SimTrackManager::SimTrackManager ( bool  iCollapsePrimaryVertices = false)

Definition at line 28 of file SimTrackManager.cc.

29  : m_trksForThisEvent(nullptr),
30  m_nVertices(0),
31  m_collapsePrimaryVertices(iCollapsePrimaryVertices),
32  lastTrack(0),
33  lastHist(0),
34  theLHCTlink(nullptr) {}
bool m_collapsePrimaryVertices
unsigned int lastHist
TrackContainer * m_trksForThisEvent
unsigned int lastTrack
const edm::LHCTransportLinkContainer * theLHCTlink

◆ ~SimTrackManager()

SimTrackManager::~SimTrackManager ( )
virtual

Definition at line 36 of file SimTrackManager.cc.

References deleteTracks(), and m_trksForThisEvent.

36  {
37  if (m_trksForThisEvent != nullptr)
38  deleteTracks();
39 }
TrackContainer * m_trksForThisEvent

◆ SimTrackManager() [2/2]

SimTrackManager::SimTrackManager ( const SimTrackManager )
delete

Member Function Documentation

◆ addTkCaloStateInfo()

void SimTrackManager::addTkCaloStateInfo ( uint32_t  t,
const std::pair< math::XYZVectorD, math::XYZTLorentzVectorD > &  p 
)
inline

Definition at line 75 of file SimTrackManager.h.

References mapTkCaloStateInfo, AlCaHLTBitMon_ParallelJobs::p, and submitPVValidationJobs::t.

Referenced by EventAction::addTkCaloStateInfo().

75  {
76  std::map<uint32_t, std::pair<math::XYZVectorD, math::XYZTLorentzVectorD> >::const_iterator it =
77  mapTkCaloStateInfo.find(t);
78 
79  if (it == mapTkCaloStateInfo.end()) {
80  mapTkCaloStateInfo.insert(std::pair<uint32_t, std::pair<math::XYZVectorD, math::XYZTLorentzVectorD> >(t, p));
81  }
82  }
std::map< uint32_t, std::pair< math::XYZVectorD, math::XYZTLorentzVectorD > > mapTkCaloStateInfo

◆ addTrack()

void SimTrackManager::addTrack ( TrackWithHistory iTrack,
bool  inHistory,
bool  withAncestor 
)
inline

Definition at line 63 of file SimTrackManager.h.

References ancestorList, idsave, m_trksForThisEvent, TrackWithHistory::parentID(), and TrackWithHistory::trackID().

Referenced by EventAction::addTrack().

63  {
64  std::pair<int, int> thePair(iTrack->trackID(), iTrack->parentID());
65  idsave.push_back(thePair);
66  if (inHistory) {
67  m_trksForThisEvent->push_back(iTrack);
68  }
69  if (withAncestor) {
70  std::pair<int, int> thisPair(iTrack->trackID(), 0);
71  ancestorList.push_back(thisPair);
72  }
73  }
int parentID() const
TrackContainer * m_trksForThisEvent
std::vector< std::pair< int, int > > idsave
std::vector< std::pair< int, int > > ancestorList
unsigned int trackID() const

◆ cleanTkCaloStateInfoMap()

void SimTrackManager::cleanTkCaloStateInfoMap ( )

Definition at line 187 of file SimTrackManager.cc.

References mapTkCaloStateInfo, and edm::swap().

Referenced by EventAction::EndOfEventAction(), and reset().

187  {
188  mapTkCaloStateInfo.clear();
189  std::map<uint32_t, std::pair<math::XYZVectorD, math::XYZTLorentzVectorD> >().swap(mapTkCaloStateInfo);
190 }
void swap(Association< C > &lhs, Association< C > &rhs)
Definition: Association.h:117
std::map< uint32_t, std::pair< math::XYZVectorD, math::XYZTLorentzVectorD > > mapTkCaloStateInfo

◆ cleanTracksWithHistory()

void SimTrackManager::cleanTracksWithHistory ( )

Definition at line 283 of file SimTrackManager.cc.

References fillMotherList(), idsave, lastTrack, LogDebug, m_trksForThisEvent, EgammaValidation_cff::num, saveTrackAndItsBranch(), submitPVValidationJobs::t, and geometryCSVtoXML::xx.

Referenced by EventAction::prepareForNewPrimary(), and storeTracks().

283  {
284  if ((*m_trksForThisEvent).empty() && idsave.empty()) {
285  return;
286  }
287 
288 #ifdef DebugLog
289  LogDebug("SimTrackManager") << "SimTrackManager::cleanTracksWithHistory has " << idsave.size()
290  << " mother-daughter relationships stored with lastTrack = " << lastTrack;
291 #endif
292 
293  if (lastTrack > 0 && lastTrack >= (*m_trksForThisEvent).size()) {
294  lastTrack = 0;
295  edm::LogError("SimTrackManager") << " SimTrackManager::cleanTracksWithHistory track index corrupted";
296  }
297 
298  stable_sort(m_trksForThisEvent->begin() + lastTrack, m_trksForThisEvent->end(), trkIDLess());
299 
300  stable_sort(idsave.begin(), idsave.end());
301 
302 #ifdef DebugLog
303  LogDebug("SimTrackManager") << " SimTrackManager::cleanTracksWithHistory knows " << m_trksForThisEvent->size()
304  << " tracks with history before branching";
305  for (unsigned int it = 0; it < (*m_trksForThisEvent).size(); it++) {
306  LogDebug("SimTrackManager") << " 1 - Track in position " << it << " G4 track number "
307  << (*m_trksForThisEvent)[it]->trackID() << " mother "
308  << (*m_trksForThisEvent)[it]->parentID() << " status "
309  << (*m_trksForThisEvent)[it]->saved();
310  }
311 #endif
312 
313  for (auto& t : *m_trksForThisEvent) {
314  if (t->saved()) {
316  }
317  }
318  unsigned int num = lastTrack;
319  for (unsigned int it = lastTrack; it < m_trksForThisEvent->size(); ++it) {
320  TrackWithHistory* t = (*m_trksForThisEvent)[it];
321  int g4ID = t->trackID();
322  if (t->saved() == true) {
323  if (it > num)
324  (*m_trksForThisEvent)[num] = t;
325  ++num;
326  for (auto& xx : idsave) {
327  if (xx.first == g4ID) {
328  xx.second = g4ID;
329  break;
330  }
331  }
332  } else {
333  delete t;
334  }
335  }
336 
337  m_trksForThisEvent->resize(num);
338 
339 #ifdef DebugLog
340  LogDebug("SimTrackManager") << " AFTER CLEANING, I GET " << (*m_trksForThisEvent).size()
341  << " tracks to be saved persistently";
342  for (unsigned int it < (*m_trksForThisEvent).size(); ++it) {
343  LogDebug("SimTrackManager") << " Track in position " << it << " G4 track number "
344  << (*m_trksForThisEvent)[it]->trackID() << " mother "
345  << (*m_trksForThisEvent)[it]->parentID() << " Status "
346  << (*m_trksForThisEvent)[it]->saved() << " id "
347  << (*m_trksForThisEvent)[it]->particleID()
348  << " E(MeV)= " << (*m_trksForThisEvent)[it]->totalEnergy();
349  }
350 #endif
351 
352  fillMotherList();
353 
354  lastTrack = (*m_trksForThisEvent).size();
355 }
Log< level::Error, false > LogError
TrackContainer * m_trksForThisEvent
std::vector< std::pair< int, int > > idsave
void saveTrackAndItsBranch(TrackWithHistory *)
this saves a track and all its parents looping over the non ordered vector
unsigned int lastTrack
#define LogDebug(id)

◆ cleanVertexMap()

void SimTrackManager::cleanVertexMap ( )
private

Definition at line 181 of file SimTrackManager.cc.

References m_nVertices, and m_vertexMap.

Referenced by reset().

181  {
182  m_vertexMap.clear();
184  m_nVertices = 0;
185 }
std::map< int, MapVertexPositionVector > MotherParticleToVertexMap
MotherParticleToVertexMap m_vertexMap

◆ deleteTracks()

void SimTrackManager::deleteTracks ( )

Definition at line 62 of file SimTrackManager.cc.

References mps_fire::i, and m_trksForThisEvent.

Referenced by EventAction::EndOfEventAction(), and ~SimTrackManager().

62  {
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 }
TrackContainer * m_trksForThisEvent

◆ fillMotherList()

void SimTrackManager::fillMotherList ( )
private

Definition at line 251 of file SimTrackManager.cc.

References ancestorList, first, idsave, idSavedTrack(), lastHist, LogDebug, and dqmiodumpmetadata::n.

Referenced by cleanTracksWithHistory().

251  {
252  if (!ancestorList.empty() && lastHist > ancestorList.size()) {
253  lastHist = ancestorList.size();
254  edm::LogError("SimTrackManager") << " SimTrackManager::fillMotherList track index corrupted";
255  }
256  /*
257  edm::LogVerbatim("SimTrackManager") << "### SimTrackManager::fillMotherList: "
258  << idsave.size() << " saved; ancestor: " << lastHist
259  << " " << ancestorList.size() << std::endl;
260  for (unsigned int i = 0; i< idsave.size(); ++i) {
261  edm::LogVerbatim("SimTrackManager") << " ISV: Track ID = " << (idsave[i]).first
262  << " Mother ID = " << (idsave[i]).second << std::endl;
263  }
264  */
265  for (unsigned int n = lastHist; n < ancestorList.size(); ++n) {
266  int theMotherId = idSavedTrack((ancestorList[n]).first);
267  ancestorList[n].second = theMotherId;
268  /*
269  std::cout << " ANC: Track ID = " << (ancestorList[n]).first
270  << " Mother ID = " << (ancestorList[n]).second << std::endl;
271  */
272 #ifdef DebugLog
273  LogDebug("SimTrackManager") << "Track ID = " << (ancestorList[n]).first
274  << " Mother ID = " << (ancestorList[n]).second;
275 #endif
276  }
277 
278  lastHist = ancestorList.size();
279 
280  idsave.clear();
281 }
Log< level::Error, false > LogError
unsigned int lastHist
std::vector< std::pair< int, int > > idsave
int idSavedTrack(int) const
std::vector< std::pair< int, int > > ancestorList
#define LogDebug(id)

◆ getOrCreateVertex()

int SimTrackManager::getOrCreateVertex ( TrackWithHistory trkH,
int  iParentID,
G4SimEvent simEvent 
)
private

Definition at line 150 of file SimTrackManager.cc.

References G4SimEvent::add(), TrackWithHistory::creatorProcess(), TrackWithHistory::globalTime(), l1ctLayer2EG_cff::id, m_nVertices, m_trksForThisEvent, m_vertexMap, class-composition::parent, OfflineHarvestingSequence_cosmic::ptype, TrackWithHistory::vertexPosition(), and geometryCSVtoXML::xx.

Referenced by reallyStoreTracks().

150  {
151  int parent = -1;
152  for (auto& trk : *m_trksForThisEvent) {
153  int id = trk->trackID();
154  if (id == iParentID) {
155  parent = id;
156  break;
157  }
158  }
159 
160  VertexMap::const_iterator iterator = m_vertexMap.find(parent);
161  if (iterator != m_vertexMap.end()) {
162  // loop over saved vertices
163  for (auto& xx : m_vertexMap[parent]) {
164  if ((trkH->vertexPosition() - xx.second).Mag2() < 1.e-6) {
165  return xx.first;
166  }
167  }
168  }
169 
170  unsigned int ptype = 0;
171  const G4VProcess* pr = trkH->creatorProcess();
172  if (nullptr != pr) {
173  ptype = pr->GetProcessSubType();
174  }
175  simEvent->add(new G4SimVertex(trkH->vertexPosition(), trkH->globalTime(), parent, ptype));
177  ++m_nVertices;
178  return (m_nVertices - 1);
179 }
double globalTime() const
const math::XYZVectorD & vertexPosition() const
TrackContainer * m_trksForThisEvent
void add(G4SimTrack *t)
Definition: G4SimEvent.h:33
MotherParticleToVertexMap m_vertexMap
std::pair< int, math::XYZVectorD > MapVertexPosition
this map contains association between vertex number and position
const G4VProcess * creatorProcess() const

◆ getTrackByID()

TrackWithHistory* SimTrackManager::getTrackByID ( unsigned int  trackID,
bool  strict = false 
) const
inline

Definition at line 104 of file SimTrackManager.h.

References Exception, m_trksForThisEvent, and HLT_2022v15_cff::track.

Referenced by CaloSD::findBoundaryCrossingParent(), and EventAction::getTrackByID().

104  {
105  bool trackFound = false;
107  if (m_trksForThisEvent == nullptr) {
108  throw cms::Exception("Unknown", "SimTrackManager") << "m_trksForThisEvent is a nullptr, cannot get any track!";
109  }
110  for (unsigned int itr = 0; itr < (*m_trksForThisEvent).size(); ++itr) {
111  if ((*m_trksForThisEvent)[itr]->trackID() == trackID) {
112  track = (*m_trksForThisEvent)[itr];
113  trackFound = true;
114  break;
115  }
116  }
117  if (!trackFound) {
118  if (strict) {
119  throw cms::Exception("Unknown", "SimTrackManager")
120  << "Attempted to get track " << trackID << " from SimTrackManager, but it's not in m_trksForThisEvent ("
121  << (*m_trksForThisEvent).size() << " tracks in m_trksForThisEvent)"
122  << "\n";
123  }
124  return nullptr;
125  }
126  return track;
127  }
TrackContainer * m_trksForThisEvent

◆ giveMotherNeeded()

int SimTrackManager::giveMotherNeeded ( int  i) const
inline

Definition at line 84 of file SimTrackManager.h.

References mps_fire::i, idsave, and edm::second().

Referenced by CaloSD::saveHit().

84  {
85  int theResult = 0;
86  for (unsigned int itr = 0; itr < idsave.size(); itr++) {
87  if ((idsave[itr]).first == i) {
88  theResult = (idsave[itr]).second;
89  break;
90  }
91  }
92  return theResult;
93  }
U second(std::pair< T, U > const &p)
std::vector< std::pair< int, int > > idsave

◆ idSavedTrack()

int SimTrackManager::idSavedTrack ( int  id) const
private

Definition at line 192 of file SimTrackManager.cc.

References ancestorList, first, l1ctLayer2EG_cff::id, globals_cff::id1, idsave, dqmiolumiharvest::j, dqmiodumpmetadata::n, notFound, and edm::second().

Referenced by fillMotherList().

192  {
193  int idMother = id;
194  if (id > 0) {
195  unsigned int n = idsave.size();
196  if (0 < n) {
197  int jmax = n - 1;
198  int j, id1;
199 
200  // first loop forward
201  bool notFound = true;
202  for (j = 0; j <= jmax; ++j) {
203  if ((idsave[j]).first == idMother) {
204  id1 = (idsave[j]).second;
205  if (0 == id1 || id1 == idMother) {
206  return id1;
207  }
208  jmax = j - 1;
209  idMother = id1;
210  notFound = false;
211  break;
212  }
213  }
214  if (notFound) {
215  return 0;
216  }
217 
218  // recursive loop
219  do {
220  notFound = true;
221  // search ID scan backward
222  for (j = jmax; j >= 0; --j) {
223  if ((idsave[j]).first == idMother) {
224  id1 = (idsave[j]).second;
225  if (0 == id1 || id1 == idMother) {
226  return id1;
227  }
228  jmax = j - 1;
229  idMother = id1;
230  notFound = false;
231  break;
232  }
233  }
234  if (notFound) {
235  // ID not in the list of saved track - look into ancestors
236  jmax = ancestorList.size() - 1;
237  for (j = jmax; j >= 0; --j) {
238  if ((ancestorList[j]).first == idMother) {
239  idMother = (ancestorList[j]).second;
240  return idMother;
241  }
242  }
243  return 0;
244  }
245  } while (!notFound);
246  }
247  }
248  return idMother;
249 }
U second(std::pair< T, U > const &p)
std::vector< std::pair< int, int > > idsave
static const GlobalPoint notFound(0, 0, 0)
std::vector< std::pair< int, int > > ancestorList

◆ operator=()

const SimTrackManager& SimTrackManager::operator= ( const SimTrackManager )
delete

◆ reallyStoreTracks()

void SimTrackManager::reallyStoreTracks ( G4SimEvent simEvent)
private

Definition at line 107 of file SimTrackManager.cc.

References G4SimEvent::add(), G4SimTrack::copyCrossedBoundaryVars(), getOrCreateVertex(), LogDebug, m_trksForThisEvent, and mapTkCaloStateInfo.

Referenced by storeTracks().

107  {
108  // loop over the (now ordered) vector and really save the tracks
109 #ifdef DebugLog
110  LogDebug("SimTrackManager") << "Inside the reallyStoreTracks method object to be stored = "
111  << m_trksForThisEvent->size();
112 #endif
113 
114  for (auto& trkH : *m_trksForThisEvent) {
115  // at this stage there is one vertex per track,
116  // so the vertex id of track N is also N
117  int ivertex = -1;
118  int ig;
119 
120  math::XYZVectorD pm(0., 0., 0.);
121  unsigned int iParentID = trkH->parentID();
122  for (auto& trk : *m_trksForThisEvent) {
123  if (trk->trackID() == iParentID) {
124  pm = trk->momentum();
125  break;
126  }
127  }
128  ig = trkH->genParticleID();
129  ivertex = getOrCreateVertex(trkH, iParentID, simEvent);
130  std::map<uint32_t, std::pair<math::XYZVectorD, math::XYZTLorentzVectorD> >::const_iterator cit =
131  mapTkCaloStateInfo.find(trkH->trackID());
132  std::pair<math::XYZVectorD, math::XYZTLorentzVectorD> tcinfo;
133  if (cit != mapTkCaloStateInfo.end()) {
134  tcinfo = cit->second;
135  }
136  G4SimTrack* g4simtrack = new G4SimTrack(trkH->trackID(),
137  trkH->particleID(),
138  trkH->momentum(),
139  trkH->totalEnergy(),
140  ivertex,
141  ig,
142  pm,
143  tcinfo.first,
144  tcinfo.second);
145  g4simtrack->copyCrossedBoundaryVars(trkH);
146  simEvent->add(g4simtrack);
147  }
148 }
void copyCrossedBoundaryVars(const TrackWithHistory *track)
Definition: G4SimTrack.h:76
int getOrCreateVertex(TrackWithHistory *, int, G4SimEvent *simEvent)
TrackContainer * m_trksForThisEvent
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::map< uint32_t, std::pair< math::XYZVectorD, math::XYZTLorentzVectorD > > mapTkCaloStateInfo
#define LogDebug(id)

◆ reset()

void SimTrackManager::reset ( void  )

Definition at line 44 of file SimTrackManager.cc.

References ancestorList, cleanTkCaloStateInfoMap(), cleanVertexMap(), mps_fire::i, idsave, lastHist, lastTrack, m_trksForThisEvent, and edm::swap().

Referenced by EventAction::BeginOfEventAction().

44  {
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 }
void cleanTkCaloStateInfoMap()
void swap(Association< C > &lhs, Association< C > &rhs)
Definition: Association.h:117
unsigned int lastHist
TrackContainer * m_trksForThisEvent
std::vector< TrackWithHistory * > TrackContainer
Definition: TrackContainer.h:8
std::vector< std::pair< int, int > > idsave
unsigned int lastTrack
std::vector< std::pair< int, int > > ancestorList

◆ resetGenID()

void SimTrackManager::resetGenID ( )
private

Definition at line 357 of file SimTrackManager.cc.

References m_trksForThisEvent, theLHCTlink, and geometryCSVtoXML::xx.

Referenced by storeTracks().

357  {
358  if (theLHCTlink == nullptr)
359  return;
360 
361  for (auto& trkH : *m_trksForThisEvent) {
362  int genParticleID = trkH->genParticleID();
363  if (genParticleID == -1) {
364  continue;
365  } else {
366  for (auto& xx : *theLHCTlink) {
367  if (xx.afterHector() == genParticleID) {
368  trkH->setGenParticleID(xx.beforeHector());
369  continue;
370  }
371  }
372  }
373  }
374 
375  theLHCTlink = nullptr;
376 }
TrackContainer * m_trksForThisEvent
const edm::LHCTransportLinkContainer * theLHCTlink

◆ saveTrackAndItsBranch()

void SimTrackManager::saveTrackAndItsBranch ( TrackWithHistory trkWHist)
private

this saves a track and all its parents looping over the non ordered vector

Definition at line 71 of file SimTrackManager.cc.

References Exception, pfDeepBoostedJetPreprocessParams_cfi::lower_bound, m_trksForThisEvent, class-composition::parent, TrackWithHistory::parentID(), and TrackWithHistory::save().

Referenced by cleanTracksWithHistory().

71  {
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  if (tk_itr != m_trksForThisEvent->end() && (*tk_itr)->trackID() == parent) {
85  saveTrackAndItsBranch(*tk_itr);
86  }
87 }
int parentID() const
Log< level::Error, false > LogError
TrackContainer * m_trksForThisEvent
void saveTrackAndItsBranch(TrackWithHistory *)
this saves a track and all its parents looping over the non ordered vector

◆ setCollapsePrimaryVertices()

void SimTrackManager::setCollapsePrimaryVertices ( bool  iSet)
inline

Definition at line 83 of file SimTrackManager.h.

References m_collapsePrimaryVertices.

Referenced by EventAction::EventAction().

83 { m_collapsePrimaryVertices = iSet; }
bool m_collapsePrimaryVertices

◆ setLHCTransportLink()

void SimTrackManager::setLHCTransportLink ( const edm::LHCTransportLinkContainer thisLHCTlink)
inline

Definition at line 128 of file SimTrackManager.h.

References theLHCTlink.

128 { theLHCTlink = thisLHCTlink; }
const edm::LHCTransportLinkContainer * theLHCTlink

◆ storeTracks()

void SimTrackManager::storeTracks ( G4SimEvent simEvent)

Definition at line 89 of file SimTrackManager.cc.

References ancestorList, cleanTracksWithHistory(), idsave, m_trksForThisEvent, reallyStoreTracks(), resetGenID(), and edm::swap().

Referenced by EventAction::EndOfEventAction().

89  {
91 
92  // fill the map with the final mother-daughter relationship
93  idsave.swap(ancestorList);
94  stable_sort(idsave.begin(), idsave.end());
95 
96  std::vector<std::pair<int, int> >().swap(ancestorList);
97 
98  // to get a backward compatible order
99  stable_sort(m_trksForThisEvent->begin(), m_trksForThisEvent->end(), trkIDLess());
100 
101  // to reset the GenParticle ID of a SimTrack to its pre-LHCTransport value
102  resetGenID();
103 
104  reallyStoreTracks(simEvent);
105 }
void swap(Association< C > &lhs, Association< C > &rhs)
Definition: Association.h:117
TrackContainer * m_trksForThisEvent
void cleanTracksWithHistory()
std::vector< std::pair< int, int > > idsave
void reallyStoreTracks(G4SimEvent *simEvent)
std::vector< std::pair< int, int > > ancestorList

◆ trackContainer()

const TrackContainer* SimTrackManager::trackContainer ( ) const
inline

Definition at line 52 of file SimTrackManager.h.

References m_trksForThisEvent.

Referenced by EventAction::trackContainer(), and CaloSD::update().

52 { return m_trksForThisEvent; }
TrackContainer * m_trksForThisEvent

◆ trackExists()

bool SimTrackManager::trackExists ( unsigned int  i) const
inline

Definition at line 94 of file SimTrackManager.h.

References RemoveAddSevLevel::flag, mps_fire::i, and m_trksForThisEvent.

Referenced by CaloSD::saveHit(), and EventAction::trackExists().

94  {
95  bool flag = false;
96  for (unsigned int itr = 0; itr < (*m_trksForThisEvent).size(); ++itr) {
97  if ((*m_trksForThisEvent)[itr]->trackID() == i) {
98  flag = true;
99  break;
100  }
101  }
102  return flag;
103  }
TrackContainer * m_trksForThisEvent

Member Data Documentation

◆ ancestorList

std::vector<std::pair<int, int> > SimTrackManager::ancestorList
private

Definition at line 154 of file SimTrackManager.h.

Referenced by addTrack(), fillMotherList(), idSavedTrack(), reset(), and storeTracks().

◆ idsave

std::vector<std::pair<int, int> > SimTrackManager::idsave
private

◆ lastHist

unsigned int SimTrackManager::lastHist
private

Definition at line 157 of file SimTrackManager.h.

Referenced by fillMotherList(), and reset().

◆ lastTrack

unsigned int SimTrackManager::lastTrack
private

Definition at line 156 of file SimTrackManager.h.

Referenced by cleanTracksWithHistory(), and reset().

◆ m_collapsePrimaryVertices

bool SimTrackManager::m_collapsePrimaryVertices
private

Definition at line 150 of file SimTrackManager.h.

Referenced by setCollapsePrimaryVertices().

◆ m_nVertices

int SimTrackManager::m_nVertices
private

Definition at line 149 of file SimTrackManager.h.

Referenced by cleanVertexMap(), and getOrCreateVertex().

◆ m_SaveSimTracks

bool SimTrackManager::m_SaveSimTracks
private

Definition at line 147 of file SimTrackManager.h.

◆ m_trksForThisEvent

TrackContainer* SimTrackManager::m_trksForThisEvent
private

◆ m_vertexMap

MotherParticleToVertexMap SimTrackManager::m_vertexMap
private

Definition at line 148 of file SimTrackManager.h.

Referenced by cleanVertexMap(), and getOrCreateVertex().

◆ mapTkCaloStateInfo

std::map<uint32_t, std::pair<math::XYZVectorD, math::XYZTLorentzVectorD> > SimTrackManager::mapTkCaloStateInfo
private

◆ theLHCTlink

const edm::LHCTransportLinkContainer* SimTrackManager::theLHCTlink
private

Definition at line 159 of file SimTrackManager.h.

Referenced by resetGenID(), and setLHCTransportLink().