9 #include "G4THitsCollection.hh"
10 #include "G4Allocator.hh"
11 #include "G4ThreeVector.hh"
13 #include "G4Colour.hh"
16 #include "gtouchable.h"
29 class GHit :
public G4VHit {
38 GHit(GTouchable* gt,
HitBitSet hbs,
const G4Step* thisStep =
nullptr,
string cScheme =
"default");
49 inline void*
operator new(size_t);
54 inline void operator delete(
void*);
69 G4Colour colour_touch, colour_hit, colour_passby;
71 bool setColorSchema();
76 GTouchable* gtouchable;
80 vector<double> edeps, times;
81 vector<G4ThreeVector> globalPositions;
82 vector<G4ThreeVector> localPositions;
87 vector<string> processNames;
90 vector<double> stepSize;
93 std::optional<double> totalEnergyDeposited;
95 G4ThreeVector avgGlobalPosition;
96 G4ThreeVector avgLocalPosition;
100 bool addHitInfosForBitIndex(
size_t bitIndex,
bool test,
const G4Step* thisStep);
108 [[nodiscard]]
inline vector<double>
getEdeps()
const {
return edeps; }
114 [[nodiscard]]
inline vector<double>
getTimes()
const {
return times; }
132 [[nodiscard]]
inline vector<int>
getPids()
const {
return pids; }
139 [[nodiscard]]
inline int getPid()
const {
return pids.front(); }
145 [[nodiscard]]
inline vector<double>
getEs()
const {
return Es; }
152 [[nodiscard]]
inline double getE()
const {
return Es.front(); }
164 [[nodiscard]]
inline GTouchable*
getGTouchable()
const {
return gtouchable; }
170 [[nodiscard]]
inline vector<GIdentifier>
getGID()
const {
return gtouchable->getIdentity(); }
210 inline void* GHit::operator
new(size_t) {
215 inline void GHit::operator
delete(
void* hit) {
Represents a hit in the detector.
GTouchable * getGTouchable() const
Returns the GTouchable associated with the hit.
void addHitInfosForBitset(HitBitSet hbs, const G4Step *thisStep)
Adds hit information based on a HitBitSet.
int getPid() const
Returns the first particle ID.
vector< G4ThreeVector > getLocalPositions() const
Gets the local positions recorded for the hit.
~GHit()=default
Destructor for GHit.
vector< int > getPids() const
Gets the particle IDs recorded for the hit.
GHit(GTouchable *gt, HitBitSet hbs, const G4Step *thisStep=nullptr, string cScheme="default")
Constructor for GHit.
vector< double > getEs() const
Gets the energy values recorded in bit 1.
vector< GIdentifier > getGID() const
Returns the detector element identity.
vector< double > getDetectorDimensions() const
Returns the dimensions of the detector element.
bool is_same_hit(GHit *hit)
Compare two ghits.
vector< double > getEdeps() const
Gets the energy depositions from each step.
string getProcessName() const
Gets the process name associated with the hit.
void Draw() override
Draws the hit using Geant4 visualization.
void calculateInfosForBit(int bit)
Calculates averaged hit information for the specified bit.
G4ThreeVector getAvgGlobaPosition()
Computes the average global position of the hit.
vector< double > getTimes() const
Gets the time stamps for each step.
double getTotalEnergyDeposited()
Computes the total energy deposited.
vector< int > getTTID()
Returns the touchable identity values as integers.
G4ThreeVector getAvgLocalPosition()
Computes the average local position of the hit.
double getE() const
Returns the first energy value.
vector< G4ThreeVector > getGlobalPositions() const
Gets the global positions recorded for the hit.
Defines the HitBitSet and associated bit indices.
std::bitset< NHITBITS > HitBitSet
G4ThreadLocal G4Allocator< GHit > * GHitAllocator
G4THitsCollection< GHit > GHitsCollection