15 #ifndef DD4HEP_GEANT4HITS_H
16 #define DD4HEP_GEANT4HITS_H
25 #include "G4StepPoint.hh"
31 namespace Simulation {
35 class Geant4TrackerHit;
36 class Geant4CalorimeterHit;
49 virtual bool operator()(
const HIT* h)
const = 0;
66 return pos == h->position;
154 Geant4TrackerHit(
int track_id,
int pdg_id,
double deposit,
double time_stamp);
165 void *
operator new(size_t);
167 void operator delete(
void *ptr);
193 void *
operator new(size_t);
195 void operator delete(
void *ptr);
201 #endif // DD4HEP_GEANT4HITS_H
int pdgID
Particle ID from the PDG table.
std::vector< MonteCarloContrib > Contributions
Position position
Hit position.
Geant4Hit()=default
Standard constructor.
Deprecated: Base class for hit comparisons.
Geant4TrackerHit & clear()
Clear hit content.
Direction momentum
Hit direction.
virtual ~Geant4CalorimeterHit()
Default destructor.
double deposit
Total energy deposit in this hit.
Deprecated: Seek the hits of an arbitrary collection for the same position.
Class of the Geant4 toolkit. See http://www-geant4.kek.jp/Reference.
double energyDeposit
Total energy deposit.
virtual ~HitCompare()
Default destructor.
Deprecated: Geant4 calorimeter hit class for deprecated sensitive detectors.
MonteCarloContrib()=default
Default constructor.
virtual ~Geant4Hit()
Default destructor.
MonteCarloContrib(int track_id, int pdg, double dep, double time_stamp)
Initializing constructor.
Geant4TrackerHit & operator=(const Geant4TrackerHit &c)
Assignment operator.
Contributions truth
Hit contributions by individual particles.
Geant4TrackerHit & storePoint(G4Step *step, G4StepPoint *point)
Store Geant4 point and step information into tracker hit structure.
static bool isGeantino(G4Track *track)
Check if the Geant4 track is a Geantino.
Contribution truth
Monte Carlo / Geant4 information.
Position position
Hit position.
ROOT::Math::XYZVector Position
HitPositionCompare(const Position &p)
Constructor.
MonteCarloContrib(int track_id, double dep, double time_stamp)
Initializing constructor.
int trackID
Geant 4 Track identifier.
Deprecated: Geant4 tracker hit class for deprecated sensitive detectors.
double length
Length of the track segment contributing to this hit.
virtual bool operator()(const HIT *h) const
Comparison function.
virtual bool operator()(const HIT *h) const =0
Comparison function.
MonteCarloContrib Contribution
virtual ~HitPositionCompare()
Default destructor.
virtual ~Geant4TrackerHit()
Default destructor.
MonteCarloContrib & operator=(const MonteCarloContrib &c)=default
Assignment operator.
double time
Timestamp when this energy was deposited.
Deprecated: basic geant4 hit class for deprecated sensitive detectors.
static Contribution extractContribution(G4Step *step)
Extract the MC contribution for a given hit from the step information.
Geant4CalorimeterHit(const Position &cell_pos)
Standard constructor.
Geant4TrackerHit()
Default constructor.