CMS 3D CMS Logo

TrackInformation.h
Go to the documentation of this file.
1 #ifndef SimG4Core_TrackInformation_H
2 #define SimG4Core_TrackInformation_H
3 
4 #include "G4VUserTrackInformation.hh"
5 #include "G4Allocator.hh"
6 #include "G4Track.hh"
9 
10 class TrackInformation : public G4VUserTrackInformation {
11 public:
13  ~TrackInformation() override = default;
14  inline void *operator new(size_t);
15  inline void operator delete(void *TrackInformation);
16 
17  bool storeTrack() const { return storeTrack_; }
19  void setStoreTrack() {
20  storeTrack_ = true;
21  isInHistory_ = true;
22  }
23 
24  bool isPrimary() const { return isPrimary_; }
25  void setPrimary(bool v) { isPrimary_ = v; }
26 
27  bool hasHits() const { return hasHits_; }
28  void setHasHits(bool v) { hasHits_ = v; }
29 
32 
33  bool isInHistory() const { return isInHistory_; }
34  void putInHistory() { isInHistory_ = true; }
35 
36  bool isAncestor() const { return flagAncestor_; }
37  void setAncestor() { flagAncestor_ = true; }
38 
39  // Calo section
40  int getIDonCaloSurface() const { return idOnCaloSurface_; }
41  void setIDonCaloSurface(int id, int ical, int last, int pdgID, double p) {
43  idCaloVolume_ = ical;
47  }
48  int getIDCaloVolume() const { return idCaloVolume_; }
49  int getIDLastVolume() const { return idLastVolume_; }
50  bool caloIDChecked() const { return caloIDChecked_; }
51  void setCaloIDChecked(bool f) { caloIDChecked_ = f; }
54  double caloSurfaceParticleP() const { return caloSurfaceParticleP_; }
56 
57  // Boundary crossing variables
58  void setCrossedBoundary(const G4Track *track);
59  bool crossedBoundary() const { return crossedBoundary_; }
62  bool startedInFineVolume() const { return startedInFineVolume_; }
63  void setStartedInFineVolume(bool flag = true) {
66  }
68 
69  // Generator information
70  int genParticlePID() const { return genParticlePID_; }
71  void setGenParticlePID(int id) { genParticlePID_ = id; }
72  double genParticleP() const { return genParticleP_; }
73  void setGenParticleP(double p) { genParticleP_ = p; }
74 
75  // remember the PID of particle entering the CASTOR detector. This is needed
76  // in order to scale the hadronic response
77  bool hasCastorHit() const { return hasCastorHit_; }
78  void setCastorHitPID(const int pid) {
79  hasCastorHit_ = true;
80  castorHitPID_ = pid;
81  }
82  int getCastorHitPID() const { return castorHitPID_; }
83 
84  // methods for MTD info management
85  //
86  void setFromTtoBTL() { mtdStatus_ |= 1 << 0; } // 1st bit
87  bool isFromTtoBTL() const { return (mtdStatus_ >> 0) & 1; }
88  void setFromBTLtoT() { mtdStatus_ |= 1 << 1; } // 2nd bit
89  bool isFromBTLtoT() const { return (mtdStatus_ >> 1) & 1; }
90  void setBTLdaughter() { mtdStatus_ |= 1 << 2; } // 3rd bit
91  bool isBTLdaughter() const { return (mtdStatus_ >> 2) & 1; }
92  void setBTLlooper() { mtdStatus_ |= 1 << 3; } // 4th bit
93  bool isBTLlooper() const { return (mtdStatus_ >> 3) & 1; }
94 
95  int idAtBTLentrance() const { return idAtBTLentrance_; }
97 
98  void Print() const override;
99 
100 private:
101  bool storeTrack_{false};
102  bool isPrimary_{false};
103  bool hasHits_{false};
105  bool isInHistory_{false};
106  bool flagAncestor_{false};
107  bool caloIDChecked_{false};
108  bool crossedBoundary_{false};
109  bool startedInFineVolume_{false};
111  bool hasCastorHit_{false};
113  int idCaloVolume_{-1};
114  int idLastVolume_{-1};
119  uint8_t mtdStatus_{0};
120  double genParticleP_{0.};
124 };
125 
126 extern G4ThreadLocal G4Allocator<TrackInformation> *fpTrackInformationAllocator;
127 
128 inline void *TrackInformation::operator new(size_t) {
130  fpTrackInformationAllocator = new G4Allocator<TrackInformation>;
131  return (void *)fpTrackInformationAllocator->MallocSingle();
132 }
133 
134 inline void TrackInformation::operator delete(void *trkInfo) {
135  fpTrackInformationAllocator->FreeSingle((TrackInformation *)trkInfo);
136 }
137 
138 #endif
bool storeTrack() const
double genParticleP() const
bool isPrimary() const
int getIDCaloVolume() const
void setStoreTrack()
can only be set to true, cannot be reset to false!
bool crossedBoundary() const
int getCastorHitPID() const
math::XYZTLorentzVectorF momentumAtBoundary_
bool hasCastorHit() const
bool isAncestor() const
const math::XYZTLorentzVectorF & getPositionAtBoundary() const
bool isBTLdaughter() const
void setGenParticlePID(int id)
G4ThreadLocal G4Allocator< TrackInformation > * fpTrackInformationAllocator
bool isInHistory() const
void setPrimary(bool v)
void setStartedInFineVolume(bool flag=true)
bool caloIDChecked() const
void setCastorHitPID(const int pid)
double caloSurfaceParticleP() const
int caloSurfaceParticlePID() const
void setCaloSurfaceParticlePID(int id)
void setIdAtBTLentrance(int id)
void setCaloSurfaceParticleP(double p)
double f[11][100]
~TrackInformation() override=default
bool hasHits() const
void setHasHits(bool v)
math::XYZTLorentzVectorF positionAtBoundary_
const math::XYZTLorentzVectorF & getMomentumAtBoundary() const
int getIDonCaloSurface() const
bool isGeneratedSecondary() const
int getIDLastVolume() const
int genParticlePID() const
bool isFromBTLtoT() const
void Print() const override
void setGenParticleP(double p)
bool isBTLlooper() const
void setGeneratedSecondary(bool v)
bool startedInFineVolumeIsSet()
bool startedInFineVolume() const
void setCaloIDChecked(bool f)
void setIDonCaloSurface(int id, int ical, int last, int pdgID, double p)
int idAtBTLentrance() const
bool isFromTtoBTL() const
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< float > > XYZTLorentzVectorF
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:22
void setCrossedBoundary(const G4Track *track)