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, HitBitSet hbs, const G4Step* thisStep = nullptr, string cScheme = "default");
39 
43  ~GHit() = default;
44 
49  inline void* operator new(size_t);
50 
54  inline void operator delete(void*);
55 
59  void Draw() override;
60 
66  bool is_same_hit(GHit* hit);
67 
68 private:
69  G4Colour colour_touch, colour_hit, colour_passby;
70 
71  bool setColorSchema();
72 
73  string colorSchema;
74 
75  // GTouchable saved here so it can be used in the overloaded == function
76  GTouchable* gtouchable;
77 
78  // hit data, selected by HitBitSet, to be collected for each step
79  // always present:
80  vector<double> edeps, times;
81  vector<G4ThreeVector> globalPositions;
82  vector<G4ThreeVector> localPositions;
83 
84  // bit 1
85  vector<int> pids;
86  vector<double> Es;
87  vector<string> processNames;
88 
89  // bit 2
90  vector<double> stepSize;
91 
92  // calculated/single quantities
93  std::optional<double> totalEnergyDeposited;
94  double averageTime;
95  G4ThreeVector avgGlobalPosition;
96  G4ThreeVector avgLocalPosition;
97  string processName;
98 
99  // build hit information based on the bit index and the touchable
100  bool addHitInfosForBitIndex(size_t bitIndex, bool test, const G4Step* thisStep);
101 
102 public:
103  // inline getters for hit information:
108  [[nodiscard]] inline vector<double> getEdeps() const { return edeps; }
109 
114  [[nodiscard]] inline vector<double> getTimes() const { return times; }
115 
120  [[nodiscard]] inline vector<G4ThreeVector> getGlobalPositions() const { return globalPositions; }
121 
126  [[nodiscard]] inline vector<G4ThreeVector> getLocalPositions() const { return localPositions; }
127 
132  [[nodiscard]] inline vector<int> getPids() const { return pids; }
133 
139  [[nodiscard]] inline int getPid() const { return pids.front(); }
140 
145  [[nodiscard]] inline vector<double> getEs() const { return Es; }
146 
152  [[nodiscard]] inline double getE() const { return Es.front(); }
153 
158  [[nodiscard]] inline string getProcessName() const { return processName; }
159 
164  [[nodiscard]] inline GTouchable* getGTouchable() const { return gtouchable; }
165 
170  [[nodiscard]] inline vector<GIdentifier> getGID() const { return gtouchable->getIdentity(); }
171 
176  [[nodiscard]] inline vector<double> getDetectorDimensions() const { return gtouchable->getDetectorDimensions(); }
177 
178  // calculated quantities
179  void calculateInfosForBit(int bit);
180 
181  double getTotalEnergyDeposited();
182 
183  double getAverageTime();
184 
185  G4ThreeVector getAvgLocalPosition();
186 
187  G4ThreeVector getAvgGlobaPosition();
188 
189  // gemc api
190  // build hit information based on the G4Step
191  void addHitInfosForBitset( HitBitSet hbs, const G4Step* thisStep);
192 
193 
201  vector<int> getTTID();
202 
203 };
204 
205 // MT definitions, as from:
206 // https://twiki.cern.ch/twiki/bin/view/Geant4/QuickMigrationGuideForGeant4V10
207 extern G4ThreadLocal G4Allocator<GHit>* GHitAllocator;
208 using GHitsCollection = G4THitsCollection<GHit>;
209 
210 inline void* GHit::operator new(size_t) {
211  if (!GHitAllocator) GHitAllocator = new G4Allocator<GHit>;
212  return (void*)GHitAllocator->MallocSingle();
213 }
214 
215 inline void GHit::operator delete(void* hit) {
216  if (!GHitAllocator) { GHitAllocator = new G4Allocator<GHit>; }
217 
218  GHitAllocator->FreeSingle((GHit*)hit);
219 }
220 
221 
222 #endif
Represents a hit in the detector.
Definition: ghit.h:29
GTouchable * getGTouchable() const
Returns the GTouchable associated with the hit.
Definition: ghit.h:164
void addHitInfosForBitset(HitBitSet hbs, const G4Step *thisStep)
Adds hit information based on a HitBitSet.
Definition: addHitInfos.cc:17
int getPid() const
Returns the first particle ID.
Definition: ghit.h:139
vector< G4ThreeVector > getLocalPositions() const
Gets the local positions recorded for the hit.
Definition: ghit.h:126
~GHit()=default
Destructor for GHit.
vector< int > getPids() const
Gets the particle IDs recorded for the hit.
Definition: ghit.h:132
GHit(GTouchable *gt, HitBitSet hbs, const G4Step *thisStep=nullptr, string cScheme="default")
Constructor for GHit.
Definition: ghit.cc:16
vector< double > getEs() const
Gets the energy values recorded in bit 1.
Definition: ghit.h:145
vector< GIdentifier > getGID() const
Returns the detector element identity.
Definition: ghit.h:170
vector< double > getDetectorDimensions() const
Returns the dimensions of the detector element.
Definition: ghit.h:176
bool is_same_hit(GHit *hit)
Compare two ghits.
Definition: ghit.cc:32
vector< double > getEdeps() const
Gets the energy depositions from each step.
Definition: ghit.h:108
string getProcessName() const
Gets the process name associated with the hit.
Definition: ghit.h:158
void Draw() override
Draws the hit using Geant4 visualization.
Definition: ghit.cc:56
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.
vector< double > getTimes() const
Gets the time stamps for each step.
Definition: ghit.h:114
double getTotalEnergyDeposited()
Computes the total energy deposited.
Definition: calculations.cc:76
double getAverageTime()
Definition: calculations.cc:87
vector< int > getTTID()
Returns the touchable identity values as integers.
Definition: ghit.cc:37
G4ThreeVector getAvgLocalPosition()
Computes the average local position of the hit.
double getE() const
Returns the first energy value.
Definition: ghit.h:152
vector< G4ThreeVector > getGlobalPositions() const
Gets the global positions recorded for the hit.
Definition: ghit.h:120
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:208