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
209
218 double averageTime;
219
225 G4ThreeVector avgGlobalPosition;
226
232 G4ThreeVector avgLocalPosition;
233
239 std::string processName;
240
252 bool addHitInfosForBitIndex(size_t bitIndex, bool test, const G4Step* thisStep);
253
257 static std::atomic<int> globalHitCounter;
258
259public:
260 // -------------------------------------------------------------------------
261 // Inline accessors (returning copies by design)
262 // -------------------------------------------------------------------------
263
268 [[nodiscard]] inline std::vector<double> getEdeps() const { return edeps; }
269
274 [[nodiscard]] inline std::vector<double> getTimes() const { return times; }
275
280 [[nodiscard]] inline std::vector<G4ThreeVector> getGlobalPositions() const { return globalPositions; }
281
286 [[nodiscard]] inline std::vector<G4ThreeVector> getLocalPositions() const { return localPositions; }
287
292 [[nodiscard]] inline std::vector<int> getPids() const { return pids; }
293
300 [[nodiscard]] inline int getPid() const { return pids.front(); }
301
306 [[nodiscard]] inline std::vector<double> getEs() const { return Es; }
307
314 [[nodiscard]] inline double getE() const { return Es.front(); }
315
323 [[nodiscard]] inline size_t nsteps() const { return Es.size(); }
324
332 [[nodiscard]] inline std::string getProcessName() const { return processName; }
333
338 [[nodiscard]] inline std::shared_ptr<GTouchable> getGTouchable() const { return gtouchable; }
339
346 [[nodiscard]] inline std::vector<GIdentifier> getGID() const { return gtouchable->getIdentity(); }
347
354 [[nodiscard]] inline std::vector<double> getDetectorDimensions() const {
355 return gtouchable->getDetectorDimensions();
356 }
357
362 [[nodiscard]] inline double getMass() const { return gtouchable->getMass(); }
363
364 // -------------------------------------------------------------------------
365 // Aggregation / calculation API
366 // -------------------------------------------------------------------------
367
378 void calculateInfosForBit(int bit);
379
387
395 double getAverageTime();
396
401 G4ThreeVector getAvgLocalPosition();
402
409 G4ThreeVector getAvgGlobaPosition();
410
411 // -------------------------------------------------------------------------
412 // Hit filling / testing helpers
413 // -------------------------------------------------------------------------
414
429 void addHitInfosForBitset(HitBitSet hbs, const G4Step* thisStep);
430
441
449 [[nodiscard]] std::vector<int> getTTID() const;
450
462 static GHit* create(const std::shared_ptr<GOptions>& gopts) {
463 HitBitSet hitBitSet;
464 auto gt = GTouchable::create(gopts);
465 auto hit = new GHit(gt, hitBitSet);
466 // Randomize between 1 and 10 steps in a deterministic, thread-safe manner.
467 hit->randomizeHitForTesting(1 + globalHitCounter.fetch_add(1, std::memory_order_relaxed) % 10);
468 return hit;
469 }
470};
471
472// MT definitions, as from:
473// https://twiki.cern.ch/twiki/bin/view/Geant4/QuickMigrationGuideForGeant4V10
474extern G4ThreadLocal G4Allocator<GHit>* GHitAllocator;
475using GHitsCollection = G4THitsCollection<GHit>;
476
477inline void* GHit::operator new(size_t) {
478 if (!GHitAllocator) GHitAllocator = new G4Allocator<GHit>;
479 return (void*)GHitAllocator->MallocSingle();
480}
481
482inline void GHit::operator delete(void* hit) {
483 if (!GHitAllocator) { GHitAllocator = new G4Allocator<GHit>; }
484
485 GHitAllocator->FreeSingle((GHit*)hit);
486}
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:292
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:300
size_t nsteps() const
Number of recorded steps for the optional-energy vector.
Definition ghit.h:323
std::vector< double > getEs() const
Get per-step total energies (when enabled).
Definition ghit.h:306
std::vector< double > getTimes() const
Get per-step global times.
Definition ghit.h:274
std::vector< GIdentifier > getGID() const
Get the detector element identity.
Definition ghit.h:346
~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:462
std::vector< G4ThreeVector > getLocalPositions() const
Get per-step local positions.
Definition ghit.h:286
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:354
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:268
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:338
double getAverageTime()
Get the average time associated with the hit.
double getMass() const
Get the sensitive element mass.
Definition ghit.h:362
G4ThreeVector getAvgLocalPosition()
Get the average local position of the hit.
double getE() const
Convenience accessor for the first energy value.
Definition ghit.h:314
std::vector< G4ThreeVector > getGlobalPositions() const
Get per-step global positions.
Definition ghit.h:280
std::string getProcessName() const
Get the representative process name for the hit.
Definition ghit.h:332
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