CMS 3D CMS Logo

TrackWithHistory.h
Go to the documentation of this file.
1 #ifndef SimG4Core_TrackWithHistory_H
2 #define SimG4Core_TrackWithHistory_H
3 
4 #include "G4Track.hh"
7 
8 #include "G4Allocator.hh"
9 
10 class G4VProcess;
17 public:
21  TrackWithHistory(const G4Track *g4track);
23 
24  inline void *operator new(size_t);
25  inline void operator delete(void *TrackWithHistory);
26 
27  void save() { saved_ = true; }
28  unsigned int trackID() const { return trackID_; }
29  int particleID() const { return particleID_; }
30  int parentID() const { return parentID_; }
31  int genParticleID() const { return genParticleID_; }
32  const math::XYZVectorD &momentum() const { return momentum_; }
33  double totalEnergy() const { return totalEnergy_; }
34  const math::XYZVectorD &vertexPosition() const { return vertexPosition_; }
35  double globalTime() const { return globalTime_; }
36  double localTime() const { return localTime_; }
37  double properTime() const { return properTime_; }
38  const G4VProcess *creatorProcess() const { return creatorProcess_; }
39  double weight() const { return weight_; }
40  void setTrackID(int i) { trackID_ = i; }
41  void setParentID(int i) { parentID_ = i; }
43  bool storeTrack() const { return storeTrack_; }
44  bool saved() const { return saved_; }
45 
46  // Boundary crossing variables
50  crossedBoundary_ = true;
51  idAtBoundary_ = id;
54  }
55  bool crossedBoundary() const { return crossedBoundary_; }
58  int getIDAtBoundary() const { return idAtBoundary_; }
64  void checkAtEnd(const G4Track *);
65 
66 private:
67  unsigned int trackID_;
69  int parentID_;
72  double totalEnergy_;
74  double globalTime_;
75  double localTime_;
76  double properTime_;
77  const G4VProcess *creatorProcess_;
78  double weight_;
80  bool saved_;
81 
82  bool isPrimary_;
87 
88  int extractGenID(const G4Track *gt) const;
89 };
90 
91 extern G4ThreadLocal G4Allocator<TrackWithHistory> *fpTrackWithHistoryAllocator;
92 
93 inline void *TrackWithHistory::operator new(size_t) {
95  fpTrackWithHistoryAllocator = new G4Allocator<TrackWithHistory>;
96  return (void *)fpTrackWithHistoryAllocator->MallocSingle();
97 }
98 
99 inline void TrackWithHistory::operator delete(void *aTwH) {
100  fpTrackWithHistoryAllocator->FreeSingle((TrackWithHistory *)aTwH);
101 }
102 
103 #endif
const G4VProcess * creatorProcess_
G4ThreadLocal G4Allocator< TrackWithHistory > * fpTrackWithHistoryAllocator
void setParentID(int i)
double globalTime() const
int parentID() const
math::XYZVectorD vertexPosition_
const math::XYZVectorD & momentum() const
math::XYZVectorD momentum_
bool storeTrack() const
double weight() const
void setCrossedBoundaryPosMom(int id, const math::XYZTLorentzVectorF position, const math::XYZTLorentzVectorF momentum)
const math::XYZVectorD & vertexPosition() const
double totalEnergy() const
void setGenParticleID(int i)
math::XYZTLorentzVectorF positionAtBoundary_
unsigned int trackID_
bool crossedBoundary() const
int genParticleID() const
int extractGenID(const G4Track *gt) const
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double > > XYZVectorD
spatial vector with cartesian internal representation
Definition: Vector3D.h:8
void setTrackID(int i)
void checkAtEnd(const G4Track *)
int getIDAtBoundary() const
math::XYZTLorentzVectorF momentumAtBoundary_
double properTime() const
bool saved() const
const math::XYZTLorentzVectorF & getMomentumAtBoundary() const
static int position[264][3]
Definition: ReadPGInfo.cc:289
double localTime() const
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< float > > XYZTLorentzVectorF
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:22
const G4VProcess * creatorProcess() const
int particleID() const
const math::XYZTLorentzVectorF & getPositionAtBoundary() const
unsigned int trackID() const
TrackWithHistory(const G4Track *g4track)