14 #ifndef DD4HEP_DDG4_GEANT4TRACKHANDLER_H
15 #define DD4HEP_DDG4_GEANT4TRACKHANDLER_H
22 #include "G4TrajectoryPoint.hh"
23 #include "G4VTouchable.hh"
24 #include "G4VSensitiveDetector.hh"
25 #include "G4ParticleDefinition.hh"
26 #include "G4DynamicParticle.hh"
27 #include "G4VProcess.hh"
32 class G4VTouchableHandle;
39 namespace Simulation {
42 class Geant4TrackHandler;
56 typedef G4ReferenceCountedHandle<G4VTouchable>
Touchable;
65 throw std::runtime_error(
"Geant4TrackHandler: NULL pointer passed to constructor!");
69 switch(
track->GetTrackStatus() ) {
70 case fAlive:
return "Alive";
71 case fStopButAlive:
return "StopButAlive";
72 case fStopAndKill:
return "StopAndKill";
73 case fKillTrackAndSecondaries:
return "KillTrackAndSecondaries";
74 case fSuspend:
return "Suspend";
75 case fPostponeToNextEvent:
return "PostponeToNextEvent";
76 default:
return "UNKNOWN";
81 operator const G4Track*()
const {
86 return track->GetTrackID();
90 return track->GetParentID();
94 return track->GetDefinition();
97 const std::string&
name()
const {
98 return trackDef()->GetParticleName();
101 const std::string&
type()
const {
102 return trackDef()->GetParticleType();
106 return track->GetPosition();
110 return track->GetVertexPosition();
114 return track->GetGlobalTime();
118 return track->GetProperTime();
122 return track->GetTotalEnergy();
126 return track->GetKineticEnergy();
130 return track->GetVelocity();
134 return track->GetTrackLength();
138 return track->GetGlobalTime();
141 G4VPhysicalVolume*
vol()
const {
142 return track->GetVolume();
145 return track->GetMomentum();
149 return track->GetNextVolume();
153 return track->GetLogicalVolumeAtVertex();
157 return track->GetTouchableHandle();
161 return track->GetNextTouchableHandle();
165 return track->GetCreatorProcess();
170 if ( p )
return p->GetProcessName();
175 return track->GetUserInformation();
178 template <
typename T>
T*
info()
const {
183 return track->GetStep();
187 return track->GetCurrentStepNumber();
191 G4ParticleDefinition* def =
trackDef();
192 return def ? def->GetPDGEncoding() : 0;
196 return track->GetDynamicParticle();
200 const G4DynamicParticle* d =
track->GetDynamicParticle();
201 if ( d )
return d->GetPrimaryParticle();
208 #endif // DD4HEP_DDG4_GEANT4TRACKHANDLER_H
G4ReferenceCountedHandle< G4VTouchable > Touchable
Geant4TrackHandler(const G4Track *t)
Initializing constructor.
const G4Step * step() const
Step information.
T * info() const
Specific user information block.
G4VPhysicalVolume * vol() const
Physical (original) volume of the track.
int parent() const
Track's parent identifier.
G4ThreeVector momentum() const
double kineticEnergy() const
Track's kinetic energy.
Helper class to ease the extraction of information from a G4Track object.
const std::string creatorName() const
Physical process of the track generation.
G4int stepNumber() const
Step number.
int pdgID() const
Access the PDG particle identification.
const Touchable & nextTouchable() const
Next touchable of the track.
const Touchable & touchable() const
Touchable of the track.
Info * userInfo() const
User information block.
const std::string & name() const
Track's particle name.
double time() const
Track time.
const G4Track * track
Reference to the track object.
const G4DynamicParticle * dynamic() const
Access the dynamic particle of the track object.
G4VUserTrackInformation Info
G4VPhysicalVolume * nextVol() const
Next physical volume of the track.
const G4ThreeVector & position() const
Track's position.
double globalTime() const
Track global time.
G4ParticleDefinition * trackDef() const
Track's particle definition.
const std::string & type() const
Track's particle type.
const char * statusName() const
Geant4TrackHandler()=delete
Inhibit default constructor.
double velocity() const
Track velocity.
double energy() const
Track's energy.
double properTime() const
Track proper time.
const G4LogicalVolume * vertexVol() const
Logical volume of the origine vertex.
const G4ThreeVector & vertex() const
Track's vertex position, where the track was created.
const G4VProcess * creatorProcess() const
Physical process of the track generation.
int id() const
Track's identifier.
const G4PrimaryParticle * primary() const
Access the primary particle of the track object (if present)
double length() const
Track length.