ghit
Loading...
Searching...
No Matches
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
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#include <atomic>
21
30class GHit : public G4VHit {
31public:
39 GHit(std::shared_ptr<GTouchable> gt, HitBitSet hbs, const G4Step* thisStep = nullptr, const std::string& cScheme = "default");
40
44 ~GHit() override = default;
45
50 inline void* operator new(size_t);
51
55 inline void operator delete(void*);
56
60 void Draw() override;
61
67 [[nodiscard]] bool is_same_hit(const GHit* hit) const;
68
69private:
70 G4Colour colour_touch, colour_hit, colour_passby;
71
72 bool setColorSchema();
73
74 std::string colorSchema;
75
76 // GTouchable saved here, so it can be used in the overloaded == function
77 std::shared_ptr<GTouchable> gtouchable;
78
79 // hit data, selected by HitBitSet, to be collected for each step
80 // always present:
81 std::vector<double> edeps, times;
82 std::vector<G4ThreeVector> globalPositions;
83 std::vector<G4ThreeVector> localPositions;
84
85 // bit 1
86 std::vector<int> pids;
87 std::vector<double> Es;
88 std::vector<std::string> processNames;
89
90 // bit 2
91 std::vector<double> stepSize;
92
93 // calculated/single quantities
94 std::optional<double> totalEnergyDeposited;
95 double averageTime;
96 G4ThreeVector avgGlobalPosition;
97 G4ThreeVector avgLocalPosition;
98 std::string processName;
99
100 // build hit information based on the bit index and the touchable
101 bool addHitInfosForBitIndex(size_t bitIndex, bool test, const G4Step* thisStep);
102
104 static std::atomic<int> globalHitCounter;
105
106public:
107 // inline getters for hit information:
112 [[nodiscard]] inline std::vector<double> getEdeps() const { return edeps; }
113
118 [[nodiscard]] inline std::vector<double> getTimes() const { return times; }
119
124 [[nodiscard]] inline std::vector<G4ThreeVector> getGlobalPositions() const { return globalPositions; }
125
130 [[nodiscard]] inline std::vector<G4ThreeVector> getLocalPositions() const { return localPositions; }
131
136 [[nodiscard]] inline std::vector<int> getPids() const { return pids; }
137
143 [[nodiscard]] inline int getPid() const { return pids.front(); }
144
149 [[nodiscard]] inline std::vector<double> getEs() const { return Es; }
150
156 [[nodiscard]] inline double getE() const { return Es.front(); }
157
158
159 [[nodiscard]] inline size_t nsteps() const { return Es.size(); }
160
161
166 [[nodiscard]] inline std::string getProcessName() const { return processName; }
167
172 [[nodiscard]] inline std::shared_ptr<GTouchable> getGTouchable() const { return gtouchable; }
173
178 [[nodiscard]] inline std::vector<GIdentifier> getGID() const { return gtouchable->getIdentity(); }
179
184 [[nodiscard]] inline std::vector<double> getDetectorDimensions() const { return gtouchable->getDetectorDimensions(); }
185
186 // calculated quantities
187 void calculateInfosForBit(int bit);
188
190
191 double getAverageTime();
192
193 G4ThreeVector getAvgLocalPosition();
194
195 G4ThreeVector getAvgGlobaPosition();
196
197 // gemc api
198 // build hit information based on the G4Step
199 void addHitInfosForBitset(HitBitSet hbs, const G4Step* thisStep);
200
202
203
211 [[nodiscard]] std::vector<int> getTTID() const;
212
213
214 // create fake GHit for testing purposes, using sector and fake dimensions
215 static GHit* create(const std::shared_ptr<GOptions>& gopts) {
216 HitBitSet hitBitSet;
217 auto gt = GTouchable::create(gopts);
218 auto hit = new GHit(gt, hitBitSet);
219 hit->randomizeHitForTesting(1 + globalHitCounter.fetch_add(1, std::memory_order_relaxed) % 10);
220 return hit;
221 }
222
223};
224
225// MT definitions, as from:
226// https://twiki.cern.ch/twiki/bin/view/Geant4/QuickMigrationGuideForGeant4V10
227extern G4ThreadLocal G4Allocator<GHit>* GHitAllocator;
228using GHitsCollection = G4THitsCollection<GHit>;
229
230inline void* GHit::operator new(size_t) {
231 if (!GHitAllocator) GHitAllocator = new G4Allocator<GHit>;
232 return (void*)GHitAllocator->MallocSingle();
233}
234
235inline void GHit::operator delete(void* hit) {
236 if (!GHitAllocator) { GHitAllocator = new G4Allocator<GHit>; }
237
238 GHitAllocator->FreeSingle((GHit*)hit);
239}
Represents a hit in the detector.
Definition ghit.h:30
GHit(std::shared_ptr< GTouchable > gt, HitBitSet hbs, const G4Step *thisStep=nullptr, const std::string &cScheme="default")
Constructor for GHit.
Definition ghit.cc:24
std::vector< int > getPids() const
Gets the particle IDs recorded for the hit.
Definition ghit.h:136
void addHitInfosForBitset(HitBitSet hbs, const G4Step *thisStep)
Adds hit information based on a HitBitSet.
int getPid() const
Returns the first particle ID.
Definition ghit.h:143
size_t nsteps() const
Definition ghit.h:159
std::vector< double > getEs() const
Gets the energy values recorded in bit 1.
Definition ghit.h:149
std::vector< double > getTimes() const
Gets the time stamps for each step.
Definition ghit.h:118
std::vector< GIdentifier > getGID() const
Returns the detector element identity.
Definition ghit.h:178
~GHit() override=default
Destructor for GHit.
static GHit * create(const std::shared_ptr< GOptions > &gopts)
Definition ghit.h:215
std::vector< G4ThreeVector > getLocalPositions() const
Gets the local positions recorded for the hit.
Definition ghit.h:130
std::vector< int > getTTID() const
Returns the touchable identity values as integers.
Definition ghit.cc:50
bool is_same_hit(const GHit *hit) const
Compare two ghits.
Definition ghit.cc:43
void randomizeHitForTesting(int nsteps)
Definition ghit.cc:115
std::vector< double > getDetectorDimensions() const
Returns the dimensions of the detector element.
Definition ghit.h:184
void Draw() override
Draws the hit using Geant4 visualization.
Definition ghit.cc:70
void calculateInfosForBit(int bit)
Calculates averaged hit information for the specified bit.
G4ThreeVector getAvgGlobaPosition()
Computes the average global position of the hit.
std::vector< double > getEdeps() const
Gets the energy depositions from each step.
Definition ghit.h:112
double getTotalEnergyDeposited()
Computes the total energy deposited.
std::shared_ptr< GTouchable > getGTouchable() const
Returns a const reference to the shared_ptr managing the GTouchable.
Definition ghit.h:172
double getAverageTime()
G4ThreeVector getAvgLocalPosition()
Computes the average local position of the hit.
double getE() const
Returns the first energy value.
Definition ghit.h:156
std::vector< G4ThreeVector > getGlobalPositions() const
Gets the global positions recorded for the hit.
Definition ghit.h:124
std::string getProcessName() const
Gets the process name associated with the hit.
Definition ghit.h:166
Defines the HitBitSet and associated bit indices.
std::bitset< NHITBITS > HitBitSet
G4ThreadLocal G4Allocator< GHit > * GHitAllocator
Definition ghit.cc:22
G4THitsCollection< GHit > GHitsCollection
Definition ghit.h:228