8#include "G4THitsCollection.hh"
9#include "G4Allocator.hh"
11#include "G4ThreeVector.hh"
58 GHit(std::shared_ptr<GTouchable> gt,
HitBitSet hbs,
const G4Step* thisStep =
nullptr,
59 const std::string& cScheme =
"default");
72 inline void*
operator new(size_t);
77 inline void operator delete(
void*);
111 G4Colour colour_touch, colour_hit, colour_passby;
118 bool setColorSchema();
121 std::string colorSchema;
144 std::vector<double> edeps;
151 std::vector<double> times;
158 std::vector<G4ThreeVector> globalPositions;
165 std::vector<G4ThreeVector> localPositions;
174 std::vector<int> pids;
181 std::vector<double> Es;
189 std::vector<std::string> processNames;
196 std::vector<double> stepSize;
207 std::optional<double> totalEnergyDeposited;
224 G4ThreeVector avgGlobalPosition;
231 G4ThreeVector avgLocalPosition;
238 std::string processName;
251 bool addHitInfosForBitIndex(
size_t bitIndex,
bool test,
const G4Step* thisStep);
256 static std::atomic<int> globalHitCounter;
267 [[nodiscard]]
inline std::vector<double>
getEdeps()
const {
return edeps; }
273 [[nodiscard]]
inline std::vector<double>
getTimes()
const {
return times; }
279 [[nodiscard]]
inline std::vector<G4ThreeVector>
getGlobalPositions()
const {
return globalPositions; }
285 [[nodiscard]]
inline std::vector<G4ThreeVector>
getLocalPositions()
const {
return localPositions; }
291 [[nodiscard]]
inline std::vector<int>
getPids()
const {
return pids; }
299 [[nodiscard]]
inline int getPid()
const {
return pids.front(); }
305 [[nodiscard]]
inline std::vector<double>
getEs()
const {
return Es; }
313 [[nodiscard]]
inline double getE()
const {
return Es.front(); }
322 [[nodiscard]]
inline size_t nsteps()
const {
return Es.size(); }
345 [[nodiscard]]
inline std::vector<GIdentifier>
getGID()
const {
return gtouchable->getIdentity(); }
442 [[nodiscard]] std::vector<int>
getTTID()
const;
455 static GHit*
create(
const std::shared_ptr<GOptions>& gopts) {
458 auto hit =
new GHit(gt, hitBitSet);
460 hit->randomizeHitForTesting(1 + globalHitCounter.fetch_add(1, std::memory_order_relaxed) % 10);
470inline void* GHit::operator
new(size_t) {
475inline void GHit::operator
delete(
void* hit) {
Stores step-by-step and aggregated information for a detector hit.
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.
std::vector< int > getPids() const
Get per-step particle PDG encodings (when enabled).
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.
size_t nsteps() const
Number of recorded steps for the optional-energy vector.
std::vector< double > getEs() const
Get per-step total energies (when enabled).
std::vector< double > getTimes() const
Get per-step global times.
std::vector< GIdentifier > getGID() const
Get the detector element identity.
~GHit() override=default
Destructor.
static GHit * create(const std::shared_ptr< GOptions > &gopts)
Create a fake hit for testing, using the current options.
std::vector< G4ThreeVector > getLocalPositions() const
Get per-step local positions.
std::vector< int > getTTID() const
Get the touchable identity values as integers.
bool is_same_hit(const GHit *hit) const
Compare this hit against another hit by sensitive-element identity.
void randomizeHitForTesting(int nsteps)
Randomize internal vectors for test-only usage.
std::vector< double > getDetectorDimensions() const
Get the sensitive-element dimensions.
void Draw() override
Visualize the hit using Geant4 visualization primitives.
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.
double getTotalEnergyDeposited()
Get the total deposited energy across all recorded steps.
std::shared_ptr< GTouchable > getGTouchable() const
Get the associated sensitive-element descriptor.
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.
std::vector< G4ThreeVector > getGlobalPositions() const
Get per-step global positions.
std::string getProcessName() const
Get the representative process name for the hit.
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
G4THitsCollection< GHit > GHitsCollection