ghit
Loading...
Searching...
No Matches
GHit Class Reference

Stores step-by-step and aggregated information for a detector hit. More...

#include <ghit.h>

Public Member Functions

 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.
 
 ~GHit () override=default
 Destructor.
 
void * operator new (size_t)
 Allocate a GHit via the G4Allocator associated with this type.
 
void operator delete (void *)
 Deallocate a GHit via the G4Allocator associated with this type.
 
void Draw () override
 Visualize the hit using Geant4 visualization primitives.
 
bool is_same_hit (const GHit *hit) const
 Compare this hit against another hit by sensitive-element identity.
 
std::vector< double > getEdeps () const
 Get per-step energy depositions.
 
std::vector< double > getTimes () const
 Get per-step global times.
 
std::vector< G4ThreeVector > getGlobalPositions () const
 Get per-step global positions.
 
std::vector< G4ThreeVector > getLocalPositions () const
 Get per-step local positions.
 
std::vector< int > getPids () const
 Get per-step particle PDG encodings (when enabled).
 
int getPid () const
 Convenience accessor for the first particle ID.
 
std::vector< double > getEs () const
 Get per-step total energies (when enabled).
 
double getE () const
 Convenience accessor for the first energy value.
 
size_t nsteps () const
 Number of recorded steps for the optional-energy vector.
 
std::string getProcessName () const
 Get the representative process name for the hit.
 
std::shared_ptr< GTouchablegetGTouchable () const
 Get the associated sensitive-element descriptor.
 
std::vector< GIdentifiergetGID () const
 Get the detector element identity.
 
std::vector< double > getDetectorDimensions () const
 Get the sensitive-element dimensions.
 
void calculateInfosForBit (int bit)
 Compute and cache derived information for the requested bit.
 
double getTotalEnergyDeposited ()
 Get the total deposited energy across all recorded steps.
 
double getAverageTime ()
 Get the average time associated with the hit.
 
G4ThreeVector getAvgLocalPosition ()
 Get the average local position of the hit.
 
G4ThreeVector getAvgGlobaPosition ()
 Get the average global position of the hit.
 
void addHitInfosForBitset (HitBitSet hbs, const G4Step *thisStep)
 Append per-step information from a G4Step according to a bitset.
 
void randomizeHitForTesting (int nsteps)
 Randomize internal vectors for test-only usage.
 
std::vector< int > getTTID () const
 Get the touchable identity values as integers.
 

Static Public Member Functions

static GHitcreate (const std::shared_ptr< GOptions > &gopts)
 Create a fake hit for testing, using the current options.
 

Detailed Description

A GHit is a G4VHit that accumulates per-step quantities while a track traverses a sensitive detector element and deposits energy.

Conceptually, this class has two layers of information:

  • Per-step vectors: always-collected quantities (energy deposition, time, local/global positions) plus optional quantities controlled by HitBitSet.
  • Aggregated quantities: totals/averages (e.g., total energy deposited, average time, average positions, representative process name) computed lazily from the per-step vectors.

The optional information is controlled by HitBitSet (see ghitConventions.h : bit meanings and expected future extensions).

Note
This class does not own the sensitive-element description. The associated GTouchable is stored as a std::shared_ptr so that the hit can be compared against other hits and can query identity/dimensions.

Definition at line 42 of file ghit.h.

Constructor & Destructor Documentation

◆ GHit()

GHit::GHit ( std::shared_ptr< GTouchable > gt,
HitBitSet hbs,
const G4Step * thisStep = nullptr,
const std::string & cScheme = "default" )

This constructor initializes the hit bookkeeping and, if thisStep is not null, immediately records per-step information for that step (both always-present data and any enabled optional data in hbs).

Parameters
gtPointer to the GTouchable describing the sensitive element producing the hit.
hbsBitset selecting which optional hit information is recorded in addition to the always-present fields.
thisStepOptional G4Step used to seed the hit with an initial step record (default: null).
cSchemeVisualization color scheme name (default: "default"). The current implementation uses a simple hard-coded scheme but keeps this field for future expansion.

Definition at line 26 of file ghit.cc.

◆ ~GHit()

GHit::~GHit ( )
overridedefault

Member Function Documentation

◆ addHitInfosForBitset()

void GHit::addHitInfosForBitset ( HitBitSet hbs,
const G4Step * thisStep )

Always records:

  • global and local positions,
  • energy deposited,
  • global time.

Then iterates over each bit of hbs and conditionally records optional information via the internal bit handler.

Parameters
hbsBitset selecting optional information to record.
thisStepStep to extract per-step values from (must not be null).

Definition at line 10 of file addHitInfos.cc.

◆ calculateInfosForBit()

void GHit::calculateInfosForBit ( int bit)

This is primarily used to compute bit-0 derived quantities (total energy, average time, average local/global positions, representative process name).

Parameters
bitBit index for which derived information should be computed.
Note
Bits beyond 0 are currently placeholders for future extensions (see ghitConventions.h : planned bits).

Implementation notes (non-Doxygen):

  • Derived quantities use a lazy-cache model (computed on first access).
  • Energy-weighted averages are used when total deposited energy is non-zero.
  • When total energy is zero, fall back to simple arithmetic averaging.

Definition at line 16 of file calculations.cc.

◆ create()

static GHit * GHit::create ( const std::shared_ptr< GOptions > & gopts)
inlinestatic

This uses GTouchable::create(gopts) to build a test touchable, constructs a hit with an empty HitBitSet, and then randomizes its contents using randomizeHitForTesting().

Parameters
goptsOptions object used to create the test GTouchable.
Returns
Newly allocated GHit pointer. Ownership is transferred to the caller.
Warning
The returned pointer must be deleted by the caller.

Definition at line 455 of file ghit.h.

◆ Draw()

void GHit::Draw ( )
override

This draws a circle at the first recorded global position and selects visual attributes based on the total energy deposited.

Note
If no visualization manager is available, or if the hit has no recorded positions, the method returns without performing any drawing.

Definition at line 65 of file ghit.cc.

◆ getAverageTime()

double GHit::getAverageTime ( )
Returns
The energy-weighted average of the time if total deposited energy is non-zero, otherwise a simple arithmetic average.
Note
The internal cache uses an "uninitialized" sentinel; computation happens on first access.

Definition at line 84 of file calculations.cc.

◆ getAvgGlobaPosition()

G4ThreeVector GHit::getAvgGlobaPosition ( )
Returns
The averaged global position energy-weighted if possible.
Note
The function name is getAvgGlobaPosition() (missing 'l' in "Global") for historical reasons.

Definition at line 105 of file calculations.cc.

◆ getAvgLocalPosition()

G4ThreeVector GHit::getAvgLocalPosition ( )
Returns
The averaged local position, energy-weighted if possible.

Definition at line 134 of file calculations.cc.

◆ getDetectorDimensions()

std::vector< double > GHit::getDetectorDimensions ( ) const
inline
Returns
A vector of doubles representing the element dimensions.
Note
This forwards to GTouchable::getDetectorDimensions().

Definition at line 353 of file ghit.h.

◆ getE()

double GHit::getE ( ) const
inline
Returns
The first energy value.
Warning
This assumes the internal Es vector is non-empty.

Definition at line 313 of file ghit.h.

◆ getEdeps()

std::vector< double > GHit::getEdeps ( ) const
inline
Returns
A copy of the vector of per-step deposited energies.

Definition at line 267 of file ghit.h.

◆ getEs()

std::vector< double > GHit::getEs ( ) const
inline
Returns
A copy of the vector of per-step energies.

Definition at line 305 of file ghit.h.

◆ getGID()

std::vector< GIdentifier > GHit::getGID ( ) const
inline
Returns
A vector of GIdentifier describing the sensitive-element identity.
Note
This forwards to GTouchable::getIdentity().

Definition at line 345 of file ghit.h.

◆ getGlobalPositions()

std::vector< G4ThreeVector > GHit::getGlobalPositions ( ) const
inline
Returns
A copy of the vector of per-step global positions.

Definition at line 279 of file ghit.h.

◆ getGTouchable()

std::shared_ptr< GTouchable > GHit::getGTouchable ( ) const
inline
Returns
A copy of the std::shared_ptr managing the GTouchable.

Definition at line 337 of file ghit.h.

◆ getLocalPositions()

std::vector< G4ThreeVector > GHit::getLocalPositions ( ) const
inline
Returns
A copy of the vector of per-step local positions.

Definition at line 285 of file ghit.h.

◆ getPid()

int GHit::getPid ( ) const
inline
Returns
The first particle ID.
Warning
This assumes the internal pids vector is non-empty.

Definition at line 299 of file ghit.h.

◆ getPids()

std::vector< int > GHit::getPids ( ) const
inline
Returns
A copy of the vector of per-step particle IDs.

Definition at line 291 of file ghit.h.

◆ getProcessName()

std::string GHit::getProcessName ( ) const
inline
Returns
The cached representative process name string.
Note
This value is typically populated when calculateInfosForBit() runs for bit 0, or when average/total getters trigger calculations.

Definition at line 331 of file ghit.h.

◆ getTimes()

std::vector< double > GHit::getTimes ( ) const
inline
Returns
A copy of the vector of per-step times.

Definition at line 273 of file ghit.h.

◆ getTotalEnergyDeposited()

double GHit::getTotalEnergyDeposited ( )
Returns
The summed energy deposition.

This method caches the result the first time it is called.

Definition at line 72 of file calculations.cc.

◆ getTTID()

vector< int > GHit::getTTID ( ) const

Converts each GIdentifier returned by getGID() into its integer value.

Returns
Vector of integer identity values (one per identifier component).

Definition at line 53 of file ghit.cc.

◆ is_same_hit()

bool GHit::is_same_hit ( const GHit * hit) const

Two hits are considered the "same" if their associated GTouchable objects compare equal (i.e. they refer to the same detector element identity according to GTouchable equality).

Parameters
hitPointer to the candidate hit to compare against (may be null).
Returns
True if hit is not null and the associated GTouchable matches.

Definition at line 46 of file ghit.cc.

◆ nsteps()

size_t GHit::nsteps ( ) const
inline
Returns
The size of the Es vector.
Note
Depending on the HitBitSet configuration, Es may remain empty even if always-present vectors have entries.

Definition at line 322 of file ghit.h.

◆ operator delete()

void GHit::operator delete ( void * hit)
inline

Definition at line 475 of file ghit.h.

◆ operator new()

void * GHit::operator new ( size_t )
inline
Returns
Pointer to the allocated storage.
Note
The allocator is thread-local (see GHitAllocator below).

Definition at line 470 of file ghit.h.

◆ randomizeHitForTesting()

void GHit::randomizeHitForTesting ( int nsteps)

Fills vectors with pseudo-random positions, times, energies, and particle IDs.

Parameters
nstepsNumber of pseudo-random steps to generate.
Warning
This is intended only for unit tests / examples. Do not use in production.

Definition at line 100 of file ghit.cc.


The documentation for this class was generated from the following files: