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
42class GHit : public G4VHit
43{
44public:
58 GHit(std::shared_ptr<GTouchable> gt, HitBitSet hbs, const G4Step* thisStep = nullptr,
59 const std::string& cScheme = "default");
60
64 ~GHit() override = default;
65
72 inline void* operator new(size_t);
73
77 inline void operator delete(void*);
78
88 void Draw() override;
89
99 [[nodiscard]] bool is_same_hit(const GHit* hit) const;
100
101private:
111 G4Colour colour_touch, colour_hit, colour_passby;
112
118 bool setColorSchema();
119
121 std::string colorSchema;
122
132 std::shared_ptr<GTouchable> gtouchable;
133
134 // -------------------------------------------------------------------------
135 // Per-step data (vectors)
136 // -------------------------------------------------------------------------
137
144 std::vector<double> edeps;
145
151 std::vector<double> times;
152
158 std::vector<G4ThreeVector> globalPositions;
159
165 std::vector<G4ThreeVector> localPositions;
166
167 // Optional per-step data, controlled by HitBitSet (see ghitConventions.h : meaning of each bit)
168
174 std::vector<int> pids;
175
181 std::vector<double> Es;
182
189 std::vector<std::string> processNames;
190
196 std::vector<double> stepSize;
197
198 // -------------------------------------------------------------------------
199 // Aggregated / calculated quantities (lazy)
200 // -------------------------------------------------------------------------
201
207 std::optional<double> totalEnergyDeposited;
208
217 double averageTime;
218
224 G4ThreeVector avgGlobalPosition;
225
231 G4ThreeVector avgLocalPosition;
232
238 std::string processName;
239
251 bool addHitInfosForBitIndex(size_t bitIndex, bool test, const G4Step* thisStep);
252
256 static std::atomic<int> globalHitCounter;
257
258public:
259 // -------------------------------------------------------------------------
260 // Inline accessors (returning copies by design)
261 // -------------------------------------------------------------------------
262
267 [[nodiscard]] inline std::vector<double> getEdeps() const { return edeps; }
268
273 [[nodiscard]] inline std::vector<double> getTimes() const { return times; }
274
279 [[nodiscard]] inline std::vector<G4ThreeVector> getGlobalPositions() const { return globalPositions; }
280
285 [[nodiscard]] inline std::vector<G4ThreeVector> getLocalPositions() const { return localPositions; }
286
291 [[nodiscard]] inline std::vector<int> getPids() const { return pids; }
292
299 [[nodiscard]] inline int getPid() const { return pids.front(); }
300
305 [[nodiscard]] inline std::vector<double> getEs() const { return Es; }
306
313 [[nodiscard]] inline double getE() const { return Es.front(); }
314
322 [[nodiscard]] inline size_t nsteps() const { return Es.size(); }
323
331 [[nodiscard]] inline std::string getProcessName() const { return processName; }
332
337 [[nodiscard]] inline std::shared_ptr<GTouchable> getGTouchable() const { return gtouchable; }
338
345 [[nodiscard]] inline std::vector<GIdentifier> getGID() const { return gtouchable->getIdentity(); }
346
353 [[nodiscard]] inline std::vector<double> getDetectorDimensions() const {
354 return gtouchable->getDetectorDimensions();
355 }
356
357 // -------------------------------------------------------------------------
358 // Aggregation / calculation API
359 // -------------------------------------------------------------------------
360
371 void calculateInfosForBit(int bit);
372
380
388 double getAverageTime();
389
394 G4ThreeVector getAvgLocalPosition();
395
402 G4ThreeVector getAvgGlobaPosition();
403
404 // -------------------------------------------------------------------------
405 // Hit filling / testing helpers
406 // -------------------------------------------------------------------------
407
422 void addHitInfosForBitset(HitBitSet hbs, const G4Step* thisStep);
423
434
442 [[nodiscard]] std::vector<int> getTTID() const;
443
455 static GHit* create(const std::shared_ptr<GOptions>& gopts) {
456 HitBitSet hitBitSet;
457 auto gt = GTouchable::create(gopts);
458 auto hit = new GHit(gt, hitBitSet);
459 // Randomize between 1 and 10 steps in a deterministic, thread-safe manner.
460 hit->randomizeHitForTesting(1 + globalHitCounter.fetch_add(1, std::memory_order_relaxed) % 10);
461 return hit;
462 }
463};
464
465// MT definitions, as from:
466// https://twiki.cern.ch/twiki/bin/view/Geant4/QuickMigrationGuideForGeant4V10
467extern G4ThreadLocal G4Allocator<GHit>* GHitAllocator;
468using GHitsCollection = G4THitsCollection<GHit>;
469
470inline void* GHit::operator new(size_t) {
471 if (!GHitAllocator) GHitAllocator = new G4Allocator<GHit>;
472 return (void*)GHitAllocator->MallocSingle();
473}
474
475inline void GHit::operator delete(void* hit) {
476 if (!GHitAllocator) { GHitAllocator = new G4Allocator<GHit>; }
477
478 GHitAllocator->FreeSingle((GHit*)hit);
479}
Stores step-by-step and aggregated information for a detector hit.
Definition ghit.h:43
GHit(std::shared_ptr< GTouchable > gt, HitBitSet hbs, const G4Step *thisStep=nullptr, const std::string &cScheme="default")
Construct a hit container and optionally seed it from a step.
Definition ghit.cc:26
std::vector< int > getPids() const
Get per-step particle PDG encodings (when enabled).
Definition ghit.h:291
void addHitInfosForBitset(HitBitSet hbs, const G4Step *thisStep)
Append per-step information from a G4Step according to a bitset.
int getPid() const
Convenience accessor for the first particle ID.
Definition ghit.h:299
size_t nsteps() const
Number of recorded steps for the optional-energy vector.
Definition ghit.h:322
std::vector< double > getEs() const
Get per-step total energies (when enabled).
Definition ghit.h:305
std::vector< double > getTimes() const
Get per-step global times.
Definition ghit.h:273
std::vector< GIdentifier > getGID() const
Get the detector element identity.
Definition ghit.h:345
~GHit() override=default
Destructor.
static GHit * create(const std::shared_ptr< GOptions > &gopts)
Create a fake hit for testing, using the current options.
Definition ghit.h:455
std::vector< G4ThreeVector > getLocalPositions() const
Get per-step local positions.
Definition ghit.h:285
std::vector< int > getTTID() const
Get the touchable identity values as integers.
Definition ghit.cc:53
bool is_same_hit(const GHit *hit) const
Compare this hit against another hit by sensitive-element identity.
Definition ghit.cc:46
void randomizeHitForTesting(int nsteps)
Randomize internal vectors for test-only usage.
Definition ghit.cc:100
std::vector< double > getDetectorDimensions() const
Get the sensitive-element dimensions.
Definition ghit.h:353
void Draw() override
Visualize the hit using Geant4 visualization primitives.
Definition ghit.cc:65
void calculateInfosForBit(int bit)
Compute and cache derived information for the requested bit.
G4ThreeVector getAvgGlobaPosition()
Get the average global position of the hit.
std::vector< double > getEdeps() const
Get per-step energy depositions.
Definition ghit.h:267
double getTotalEnergyDeposited()
Get the total deposited energy across all recorded steps.
std::shared_ptr< GTouchable > getGTouchable() const
Get the associated sensitive-element descriptor.
Definition ghit.h:337
double getAverageTime()
Get the average time associated with the hit.
G4ThreeVector getAvgLocalPosition()
Get the average local position of the hit.
double getE() const
Convenience accessor for the first energy value.
Definition ghit.h:313
std::vector< G4ThreeVector > getGlobalPositions() const
Get per-step global positions.
Definition ghit.h:279
std::string getProcessName() const
Get the representative process name for the hit.
Definition ghit.h:331
static std::shared_ptr< GTouchable > create(const std::shared_ptr< GOptions > &gopt)
Defines the hit information selection bitset for GHit.
std::bitset< NHITBITS > HitBitSet
Bitset selecting which optional hit information is recorded.
G4ThreadLocal G4Allocator< GHit > * GHitAllocator
Definition ghit.cc:22
G4THitsCollection< GHit > GHitsCollection