ghit
ghit.h
Go to the documentation of this file.
1 #pragma once
2 
3 // HitBitSet definition
4 #include "ghitConventions.h"
5 
6 // geant4
7 #include "G4VHit.hh"
8 #include "G4THitsCollection.hh"
9 #include "G4Allocator.hh"
10 #include "G4ThreeVector.hh"
11 #include "G4Step.hh"
12 #include "G4Colour.hh"
13 
14 // gemc
15 #include "gtouchable.h"
16 
17 // c++
18 #include <optional>
19 
28 class GHit : public G4VHit {
29 public:
37  GHit(GTouchable* gt, HitBitSet hbs, const G4Step* thisStep = nullptr, std::string cScheme = "default");
38 
42  ~GHit() override = default;
43 
48  inline void* operator new(size_t);
49 
53  inline void operator delete(void*);
54 
58  void Draw() override;
59 
65  bool is_same_hit(GHit* hit);
66 
67 private:
68  G4Colour colour_touch, colour_hit, colour_passby;
69 
70  bool setColorSchema();
71 
72  std::string colorSchema;
73 
74  // GTouchable saved here so it can be used in the overloaded == function
75  GTouchable* gtouchable;
76 
77  // hit data, selected by HitBitSet, to be collected for each step
78  // always present:
79  std::vector<double> edeps, times;
80  std::vector<G4ThreeVector> globalPositions;
81  std::vector<G4ThreeVector> localPositions;
82 
83  // bit 1
84  std::vector<int> pids;
85  std::vector<double> Es;
86  std::vector<std::string> processNames;
87 
88  // bit 2
89  std::vector<double> stepSize;
90 
91  // calculated/single quantities
92  std::optional<double> totalEnergyDeposited;
93  double averageTime;
94  G4ThreeVector avgGlobalPosition;
95  G4ThreeVector avgLocalPosition;
96  std::string processName;
97 
98  // build hit information based on the bit index and the touchable
99  bool addHitInfosForBitIndex(size_t bitIndex, bool test, const G4Step* thisStep);
100 
101 
102 public:
103  // inline getters for hit information:
108  [[nodiscard]] inline std::vector<double> getEdeps() const { return edeps; }
109 
114  [[nodiscard]] inline std::vector<double> getTimes() const { return times; }
115 
120  [[nodiscard]] inline std::vector<G4ThreeVector> getGlobalPositions() const { return globalPositions; }
121 
126  [[nodiscard]] inline std::vector<G4ThreeVector> getLocalPositions() const { return localPositions; }
127 
132  [[nodiscard]] inline std::vector<int> getPids() const { return pids; }
133 
139  [[nodiscard]] inline int getPid() const { return pids.front(); }
140 
145  [[nodiscard]] inline std::vector<double> getEs() const { return Es; }
146 
152  [[nodiscard]] inline double getE() const { return Es.front(); }
153 
158  [[nodiscard]] inline std::string getProcessName() const { return processName; }
159 
164  [[nodiscard]] inline GTouchable* getGTouchable() const { return gtouchable; }
165 
170  [[nodiscard]] inline std::vector<GIdentifier> getGID() const { return gtouchable->getIdentity(); }
171 
176  [[nodiscard]] inline std::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  void randomizeHitForTesting(int nsteps);
194 
195 
203  std::vector<int> getTTID();
204 
205 };
206 
207 // MT definitions, as from:
208 // https://twiki.cern.ch/twiki/bin/view/Geant4/QuickMigrationGuideForGeant4V10
209 extern G4ThreadLocal G4Allocator<GHit>* GHitAllocator;
210 using GHitsCollection = G4THitsCollection<GHit>;
211 
212 inline void* GHit::operator new(size_t) {
213  if (!GHitAllocator) GHitAllocator = new G4Allocator<GHit>;
214  return (void*)GHitAllocator->MallocSingle();
215 }
216 
217 inline void GHit::operator delete(void* hit) {
218  if (!GHitAllocator) { GHitAllocator = new G4Allocator<GHit>; }
219 
220  GHitAllocator->FreeSingle((GHit*)hit);
221 }
Represents a hit in the detector.
Definition: ghit.h:28
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
GHit(GTouchable *gt, HitBitSet hbs, const G4Step *thisStep=nullptr, std::string cScheme="default")
Constructor for GHit.
Definition: ghit.cc:37
int getPid() const
Returns the first particle ID.
Definition: ghit.h:139
std::vector< double > getDetectorDimensions() const
Returns the dimensions of the detector element.
Definition: ghit.h:176
std::vector< double > getEdeps() const
Gets the energy depositions from each step.
Definition: ghit.h:108
std::vector< GIdentifier > getGID() const
Returns the detector element identity.
Definition: ghit.h:170
std::vector< double > getTimes() const
Gets the time stamps for each step.
Definition: ghit.h:114
~GHit() override=default
Destructor for GHit.
void randomizeHitForTesting(int nsteps)
Definition: ghit.cc:22
std::vector< G4ThreeVector > getGlobalPositions() const
Gets the global positions recorded for the hit.
Definition: ghit.h:120
bool is_same_hit(GHit *hit)
Compare two ghits.
Definition: ghit.cc:53
void Draw() override
Draws the hit using Geant4 visualization.
Definition: ghit.cc:78
std::vector< double > getEs() const
Gets the energy values recorded in bit 1.
Definition: ghit.h:145
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.
double getTotalEnergyDeposited()
Computes the total energy deposited.
Definition: calculations.cc:76
double getAverageTime()
Definition: calculations.cc:87
std::vector< int > getTTID()
Returns the touchable identity values as integers.
Definition: ghit.cc:58
G4ThreeVector getAvgLocalPosition()
Computes the average local position of the hit.
double getE() const
Returns the first energy value.
Definition: ghit.h:152
std::vector< int > getPids() const
Gets the particle IDs recorded for the hit.
Definition: ghit.h:132
std::string getProcessName() const
Gets the process name associated with the hit.
Definition: ghit.h:158
std::vector< G4ThreeVector > getLocalPositions() const
Gets the local positions recorded for the hit.
Definition: ghit.h:126
Defines the HitBitSet and associated bit indices.
std::bitset< NHITBITS > HitBitSet
G4ThreadLocal G4Allocator< GHit > * GHitAllocator
Definition: ghit.cc:20
G4THitsCollection< GHit > GHitsCollection
Definition: ghit.h:210