|
DD4hep
1.30.0
Detector Description Toolkit for High Energy Physics
|
Go to the documentation of this file.
13 #ifndef DDG4_GEANT4TRACKHANDLER_H
14 #define DDG4_GEANT4TRACKHANDLER_H
21 #include <G4TrajectoryPoint.hh>
22 #include <G4VTouchable.hh>
23 #include <G4VSensitiveDetector.hh>
24 #include <G4ParticleDefinition.hh>
25 #include <G4DynamicParticle.hh>
26 #include <G4VProcess.hh>
32 class G4VTouchableHandle;
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();
142 G4ParticleDefinition* def =
trackDef();
143 return def ? def->GetPDGCharge() : 0;
147 G4ParticleDefinition* def =
trackDef();
148 return def ? def->GetPDGMass() : 0;
151 G4VPhysicalVolume*
vol()
const {
152 return track->GetVolume();
155 return track->GetMomentum();
159 return track->GetNextVolume();
163 return track->GetLogicalVolumeAtVertex();
167 return track->GetTouchableHandle();
171 return track->GetNextTouchableHandle();
175 return track->GetCreatorProcess();
180 if ( p )
return p->GetProcessName();
185 return track->GetUserInformation();
188 template <
typename T> T*
info()
const {
193 return track->GetStep();
197 return track->GetCurrentStepNumber();
201 G4ParticleDefinition* def =
trackDef();
202 return def ? def->GetPDGEncoding() : 0;
206 return track->GetDynamicParticle();
210 const G4DynamicParticle* d =
track->GetDynamicParticle();
211 if ( d )
return d->GetPrimaryParticle();
218 #endif // DDG4_GEANT4TRACKHANDLER_H
Geant4TrackHandler()=delete
Inhibit default constructor.
G4int stepNumber() const
Step number.
G4VPhysicalVolume * vol() const
Physical (original) volume of the track.
const Touchable & touchable() const
Touchable of the track.
double kineticEnergy() const
Track's kinetic energy.
const G4ThreeVector & position() const
Track's position.
double globalTime() const
Track global time.
const G4Step * step() const
Step information.
Geant4TrackHandler(const G4Track *t)
Initializing constructor.
int parent() const
Track's parent identifier.
double time() const
Track time.
G4VUserTrackInformation Info
const G4ThreeVector & vertex() const
Track's vertex position, where the track was created.
const std::string & name() const
Track's particle name.
double charge() const
Track charge.
const G4VProcess * creatorProcess() const
Physical process of the track generation.
Helper class to ease the extraction of information from a G4Track object.
const G4Track * track
Reference to the track object.
int pdgID() const
Access the PDG particle identification.
double velocity() const
Track velocity.
T * info() const
Specific user information block.
const char * statusName() const
G4ThreeVector momentum() const
double length() const
Track length.
Namespace for the AIDA detector description toolkit.
G4ParticleDefinition * trackDef() const
Track's particle definition.
const G4LogicalVolume * vertexVol() const
Logical volume of the origine vertex.
Info * userInfo() const
User information block.
const Touchable & nextTouchable() const
Next touchable of the track.
const G4PrimaryParticle * primary() const
Access the primary particle of the track object (if present)
double properTime() const
Track proper time.
const G4DynamicParticle * dynamic() const
Access the dynamic particle of the track object.
double mass() const
Track charge.
int id() const
Track's identifier.
G4ReferenceCountedHandle< G4VTouchable > Touchable
double energy() const
Track's energy.
const std::string & type() const
Track's particle type.
const std::string creatorName() const
Physical process of the track generation.
G4VPhysicalVolume * nextVol() const
Next physical volume of the track.