ghit
ghit.h
Go to the documentation of this file.
1 #ifndef GHIT_H
2 #define GHIT_H 1
3 
4 // HitBitSet definition
5 #include "ghitConventions.h"
6 
7 // geant4
8 #include "G4VHit.hh"
9 #include "G4THitsCollection.hh"
10 #include "G4Allocator.hh"
11 #include "G4ThreeVector.hh"
12 #include "G4Step.hh"
13 #include "G4Colour.hh"
14 
15 // gemc
16 #include "gtouchable.h"
17 
18 // c++
19 #include <optional>
20 
29 class GHit : public G4VHit {
30 public:
38  GHit(GTouchable *gt, const HitBitSet hbs, const G4Step *thisStep = nullptr, const string cScheme = "default");
39 
43  ~GHit() = default;
44 
50  inline void *operator new(size_t);
51 
56  inline void operator delete(void *);
57 
61  void Draw() override;
62 
68  bool is_same_hit(GHit *hit);
69 
70 private:
71 
72  G4Colour colour_touch, colour_hit, colour_passby;
73 
74  bool setColorSchema();
75 
76  string colorSchema;
77 
78  // GTouchable saved here so it can be used in the overloaded == function
79  GTouchable *gtouchable;
80 
81  // hit data, selected by HitBitSet, to be collected for each step
82  // always present:
83  vector<float> edeps, times;
84  vector <G4ThreeVector> globalPositions;
85  vector <G4ThreeVector> localPositions;
86 
87  // bit 1
88  vector<int> pids;
89  vector<float> Es;
90  vector <string> processNames;
91 
92  // bit 2
93  vector<float> stepSize;
94 
95  // calculated/single quantities
96  std::optional<float> totalEnergyDeposited;
97  float averageTime;
98  G4ThreeVector avgGlobalPosition;
99  G4ThreeVector avgLocalPosition;
100  string processName;
101 
102  // build hit information based on the bit index and the touchable
103  bool addHitInfosForBitIndex(size_t bitIndex, const bool test, const G4Step *thisStep);
104 
105 public:
106  // inline getters for hit information:
111  inline const vector<float> getEdeps() const { return edeps; }
112 
117  inline const vector<float> getTimes() const { return times; }
118 
123  inline const vector <G4ThreeVector> getGlobalPositions() const { return globalPositions; }
124 
129  inline const vector <G4ThreeVector> getLocalPositions() const { return localPositions; }
130 
135  inline const vector<int> getPids() const { return pids; }
136 
142  inline int getPid() const { return pids.front(); }
143 
148  inline vector<float> getEs() const { return Es; }
149 
155  inline float getE() const { return Es.front(); }
156 
161  inline const string getProcessName() const { return processName; }
162 
167  inline const GTouchable *getGTouchable() const { return gtouchable; }
168 
173  inline const vector <GIdentifier> getGID() const { return gtouchable->getIdentity(); }
174 
179  inline const vector<double> getDetectorDimensions() const { return gtouchable->getDetectorDimensions(); }
180 
181  // calculated quantities
182  void calculateInfosForBit(int bit);
183 
184  float getTotalEnergyDeposited();
185 
186  float getAverageTime();
187 
188  G4ThreeVector getAvgLocalPosition();
189 
190  G4ThreeVector getAvgGlobaPosition();
191 
192  // gemc api
193  // build hit information based on the G4Step
194  void addHitInfosForBitset(const HitBitSet hbs, const G4Step *thisStep);
195 
196 
204  vector<int> getTTID();
205 
206 };
207 
208 // MT definitions, as from:
209 // https://twiki.cern.ch/twiki/bin/view/Geant4/QuickMigrationGuideForGeant4V10
210 extern G4ThreadLocal G4Allocator<GHit>* GHitAllocator;
211 using GHitsCollection = G4THitsCollection<GHit>;
212 
213 inline void *GHit::operator new(size_t) {
214  if (!GHitAllocator) GHitAllocator = new G4Allocator<GHit>;
215  return (void *) GHitAllocator->MallocSingle();
216 
217 }
218 
219 inline void GHit::operator delete(void *hit) {
220  if (!GHitAllocator) {
221  GHitAllocator = new G4Allocator<GHit>;
222  }
223 
224  GHitAllocator->FreeSingle((GHit *) hit);
225 }
226 
227 
228 #endif
Represents a hit in the detector.
Definition: ghit.h:29
float getAverageTime()
Definition: calculations.cc:87
int getPid() const
Returns the first particle ID.
Definition: ghit.h:142
const vector< double > getDetectorDimensions() const
Returns the dimensions of the detector element.
Definition: ghit.h:179
const vector< G4ThreeVector > getGlobalPositions() const
Gets the global positions recorded for the hit.
Definition: ghit.h:123
~GHit()=default
Destructor for GHit.
const vector< GIdentifier > getGID() const
Returns the detector element identity.
Definition: ghit.h:173
const vector< float > getEdeps() const
Gets the energy depositions from each step.
Definition: ghit.h:111
const vector< float > getTimes() const
Gets the time stamps for each step.
Definition: ghit.h:117
const string getProcessName() const
Gets the process name associated with the hit.
Definition: ghit.h:161
GHit(GTouchable *gt, const HitBitSet hbs, const G4Step *thisStep=nullptr, const string cScheme="default")
Constructor for GHit.
Definition: ghit.cc:16
void addHitInfosForBitset(const HitBitSet hbs, const G4Step *thisStep)
Adds hit information based on a HitBitSet.
Definition: addHitInfos.cc:17
bool is_same_hit(GHit *hit)
Compare two ghits.
Definition: ghit.cc:32
void Draw() override
Draws the hit using Geant4 visualization.
Definition: ghit.cc:56
const vector< int > getPids() const
Gets the particle IDs recorded for the hit.
Definition: ghit.h:135
const vector< G4ThreeVector > getLocalPositions() const
Gets the local positions recorded for the hit.
Definition: ghit.h:129
void calculateInfosForBit(int bit)
Calculates averaged hit information for the specified bit.
Definition: calculations.cc:17
G4ThreeVector getAvgGlobaPosition()
Computes the average global position of the hit.
float getE() const
Returns the first energy value.
Definition: ghit.h:155
float getTotalEnergyDeposited()
Computes the total energy deposited.
Definition: calculations.cc:76
vector< int > getTTID()
Returns the touchable identity values as integers.
Definition: ghit.cc:37
G4ThreeVector getAvgLocalPosition()
Computes the average local position of the hit.
vector< float > getEs() const
Gets the energy values recorded in bit 1.
Definition: ghit.h:148
const GTouchable * getGTouchable() const
Returns the GTouchable associated with the hit.
Definition: ghit.h:167
Defines the HitBitSet and associated bit indices.
std::bitset< NHITBITS > HitBitSet
G4ThreadLocal G4Allocator< GHit > * GHitAllocator
Definition: ghit.cc:14
G4THitsCollection< GHit > GHitsCollection
Definition: ghit.h:211