EDepSimTrajectory.hh
Go to the documentation of this file.
1 #ifndef EDepSim_Trajectory_hh_seen
2 #define EDepSim_Trajectory_hh_seen
3 
4 #include <stdlib.h>
5 #include <vector>
6 
7 #include <globals.hh>
8 #include <G4VTrajectory.hh>
9 #include <G4Allocator.hh>
10 #include <G4ios.hh>
11 #include <G4ParticleDefinition.hh>
12 #include <G4Track.hh>
13 #include <G4Step.hh>
14 
16 
17 typedef std::vector<G4VTrajectoryPoint*> TrajectoryPointContainer;
18 
19 namespace EDepSim {class Trajectory;}
20 
21 /// EDepSim specific trajectory class to save internal bookkeeping
22 /// information. This keeps track of information about a particular particle
23 /// track. Important information kept includes how much energy this
24 /// trajectory has deposited in sensitive detectors, and it can query it's
25 /// children to find out how much energy they deposited. It also keeps track
26 /// of the process that created the trajectory, the total length in the
27 /// sensitive detectors, and points along the trajectory.
28 class EDepSim::Trajectory : public G4VTrajectory {
29 public:
30  Trajectory();
31 
32  Trajectory(const G4Track* aTrack);
34  virtual ~Trajectory();
35 
36  inline void* operator new(size_t);
37  inline void operator delete(void*);
38  inline int operator == (const EDepSim::Trajectory& right) const {
39  return (this==&right);
40  }
41 
42  /// Get the track id described by this trajectory.
43  inline G4int GetTrackID() const {return fTrackID;}
44 
45  /// Get the track id of the parent of the track described by this
46  /// trajectory.
47  inline G4int GetParentID() const {return fParentID;}
48 
49  /// Get the particle name.
50  inline G4String GetParticleName() const {return fParticleName;}
51 
52  /// Get the particle charge.
53  inline G4double GetCharge() const {return fPDGCharge;}
54 
55  /// Get the PDG MC particle number for this particle.
56  inline G4int GetPDGEncoding() const {return fPDGEncoding;}
57 
58  /// Get the interaction process that created the trajectory.
59  inline G4String GetProcessName() const {return fProcessName;}
60 
61  /// Get the initial momentum of the particle that created this trajectory.
62  inline G4ThreeVector GetInitialMomentum() const {return fInitialMomentum;}
63 
64  /// Get the initial kinetic energy of the particle that created this
65  /// trajectory.
66  G4double GetInitialKineticEnergy() const;
67 
68  /// Get the amount of energy deposited into a sensitive detector.
69  G4double GetSDEnergyDeposit() const {return fSDEnergyDeposit;}
70 
71  /// Add energy deposited into a sensitive detector.
72  void AddSDEnergyDeposit(double energy) {
75  }
76 
77  /// Get the total length of this trajectory that is in a sensitive detector.
78  G4double GetSDLength() const {return fSDLength;}
79 
80  /// Add the length that has been deposited into a sensitive detector.
81  void AddSDLength(double len) {
82  fSDLength += len;
83  }
84 
85  /// Get the total amount of energy deposited into a sensitive detector by
86  /// this trajectory and all of it's daughters.
88 
89  /// Add energy deposited into a sensitive detector by a daughter.
92  }
93 
94  /// Get the range of the particle that created this trajectory.
95  G4double GetRange() const;
96 
97  /// Mark this trajectory as one that should be saved in the output.
98  void MarkTrajectory(bool save=true);
99 
100  /// Check if this trajectory should be saved.
101  bool SaveTrajectory() const { return fSaveTrajectory;}
102 
103  // Add a new step to the trajectory.
104  virtual void AppendStep(const G4Step* aStep);
105 
106  /// Get the number of trajectory points saved with this trajectory.
107  virtual int GetPointEntries() const {return fPositionRecord->size();}
108 
109  /// Get a particular trajectory point.
110  virtual G4VTrajectoryPoint* GetPoint(G4int i) const {
111  return (*fPositionRecord)[i];
112  }
113 
114  virtual void MergeTrajectory(G4VTrajectory* secondTrajectory);
115 
116  /// Get the full definition of the particle.
117  G4ParticleDefinition* GetParticleDefinition() const;
118 
119  virtual const std::map<G4String,G4AttDef>* GetAttDefs() const;
120  virtual std::vector<G4AttValue>* CreateAttValues() const;
121 
122 private:
123 
125  G4int fTrackID;
126  G4int fParentID;
128  G4double fPDGCharge;
129  G4String fParticleName;
130  G4String fProcessName;
131  G4ThreeVector fInitialMomentum;
134  G4double fSDLength;
136 };
137 
138 #if defined G4TRACKING_ALLOC_EXPORT
139 extern G4DLLEXPORT G4Allocator<EDepSim::Trajectory> aTrajAllocator;
140 #else
141 extern G4DLLIMPORT G4Allocator<EDepSim::Trajectory> aTrajAllocator;
142 #endif
143 
144 inline void* EDepSim::Trajectory::operator new(size_t) {
145  void* aTrajectory = (void*) aTrajAllocator.MallocSingle();
146  return aTrajectory;
147 }
148 
149 inline void EDepSim::Trajectory::operator delete(void* aTrajectory) {
150  aTrajAllocator.FreeSingle((EDepSim::Trajectory*)aTrajectory);
151 }
152 #endif
G4double GetInitialKineticEnergy() const
G4ParticleDefinition * GetParticleDefinition() const
Get the full definition of the particle.
TrajectoryPointContainer * fPositionRecord
G4int GetTrackID() const
Get the track id described by this trajectory.
G4double GetRange() const
Get the range of the particle that created this trajectory.
G4DLLIMPORT G4Allocator< EDepSim::Trajectory > aTrajAllocator
G4double GetSDEnergyDeposit() const
Get the amount of energy deposited into a sensitive detector.
virtual void AppendStep(const G4Step *aStep)
G4double GetCharge() const
Get the particle charge.
G4ThreeVector GetInitialMomentum() const
Get the initial momentum of the particle that created this trajectory.
G4double GetSDLength() const
Get the total length of this trajectory that is in a sensitive detector.
virtual void MergeTrajectory(G4VTrajectory *secondTrajectory)
G4ThreeVector fInitialMomentum
virtual G4VTrajectoryPoint * GetPoint(G4int i) const
Get a particular trajectory point.
void MarkTrajectory(bool save=true)
Mark this trajectory as one that should be saved in the output.
G4String GetParticleName() const
Get the particle name.
G4String GetProcessName() const
Get the interaction process that created the trajectory.
void AddSDLength(double len)
Add the length that has been deposited into a sensitive detector.
Construct a module from components.
Definition: TG4HitSegment.h:10
bool SaveTrajectory() const
Check if this trajectory should be saved.
G4int GetPDGEncoding() const
Get the PDG MC particle number for this particle.
G4double GetSDTotalEnergyDeposit() const
def save(obj, fname)
Definition: root.py:19
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
std::vector< G4VTrajectoryPoint * > TrajectoryPointContainer
virtual std::vector< G4AttValue > * CreateAttValues() const
void AddSDEnergyDeposit(double energy)
Add energy deposited into a sensitive detector.
virtual int GetPointEntries() const
Get the number of trajectory points saved with this trajectory.
int operator==(const EDepSim::Trajectory &right) const
G4int GetParentID() const
void AddSDDaughterEnergyDeposit(double energy)
Add energy deposited into a sensitive detector by a daughter.