DD4hep - The AIDA detector description toolkit for high energy physics experiments
DD4hep  Rev:Unversioneddirectory
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Geant4TrackHandler.h
Go to the documentation of this file.
1 // $Id: $
2 //==========================================================================
3 // AIDA Detector description implementation for LCD
4 //--------------------------------------------------------------------------
5 // Copyright (C) Organisation europeenne pour la Recherche nucleaire (CERN)
6 // All rights reserved.
7 //
8 // For the licensing terms see $DD4hepINSTALL/LICENSE.
9 // For the list of contributors see $DD4hepINSTALL/doc/CREDITS.
10 //
11 // Author : M.Frank
12 //
13 //==========================================================================
14 #ifndef DD4HEP_DDG4_GEANT4TRACKHANDLER_H
15 #define DD4HEP_DDG4_GEANT4TRACKHANDLER_H
16 
17 // Framework include files
18 #include "DDG4/Defs.h"
19 
20 // Geant4 include files
21 #include "G4Track.hh"
22 #include "G4TrajectoryPoint.hh"
23 #include "G4VTouchable.hh"
24 #include "G4VSensitiveDetector.hh"
25 #include "G4ParticleDefinition.hh"
26 #include "G4DynamicParticle.hh"
27 #include "G4VProcess.hh"
28 
29 #include <stdexcept>
30 
31 // Forward declarations
32 class G4VTouchableHandle;
34 
36 namespace DD4hep {
37 
39  namespace Simulation {
40 
41  // Forward declarations;
42  class Geant4TrackHandler;
43 
45 
54  public:
56  typedef G4ReferenceCountedHandle<G4VTouchable> Touchable;
58  const G4Track* track;
60  Geant4TrackHandler() = delete;
62  Geant4TrackHandler(const G4Track* t) : track(t) {
64  if ( 0 == t ) {
65  throw std::runtime_error("Geant4TrackHandler: NULL pointer passed to constructor!");
66  }
67  }
68  const char* statusName() const {
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";
77  }
78  }
79 
81  operator const G4Track*() const {
82  return track;
83  }
85  int id() const {
86  return track->GetTrackID();
87  }
89  int parent() const {
90  return track->GetParentID();
91  }
93  G4ParticleDefinition* trackDef() const {
94  return track->GetDefinition();
95  }
97  const std::string& name() const {
98  return trackDef()->GetParticleName();
99  }
101  const std::string& type() const {
102  return trackDef()->GetParticleType();
103  }
105  const G4ThreeVector& position() const {
106  return track->GetPosition();
107  }
109  const G4ThreeVector& vertex() const {
110  return track->GetVertexPosition();
111  }
113  double globalTime() const {
114  return track->GetGlobalTime();
115  }
117  double properTime() const {
118  return track->GetProperTime();
119  }
121  double energy() const {
122  return track->GetTotalEnergy();
123  }
125  double kineticEnergy() const {
126  return track->GetKineticEnergy();
127  }
129  double velocity() const {
130  return track->GetVelocity();
131  }
133  double length() const {
134  return track->GetTrackLength();
135  }
137  double time() const {
138  return track->GetGlobalTime();
139  }
141  G4VPhysicalVolume* vol() const {
142  return track->GetVolume();
143  }
144  G4ThreeVector momentum() const {
145  return track->GetMomentum();
146  }
148  G4VPhysicalVolume* nextVol() const {
149  return track->GetNextVolume();
150  }
152  const G4LogicalVolume* vertexVol() const {
153  return track->GetLogicalVolumeAtVertex();
154  }
156  const Touchable& touchable() const {
157  return track->GetTouchableHandle();
158  }
160  const Touchable& nextTouchable() const {
161  return track->GetNextTouchableHandle();
162  }
164  const G4VProcess* creatorProcess() const {
165  return track->GetCreatorProcess();
166  }
168  const std::string creatorName() const {
169  const G4VProcess* p = creatorProcess();
170  if ( p ) return p->GetProcessName();
171  return "";
172  }
174  Info* userInfo() const {
175  return track->GetUserInformation();
176  }
178  template <typename T> T* info() const {
179  return (T*) userInfo();
180  }
182  const G4Step* step() const {
183  return track->GetStep();
184  }
186  G4int stepNumber() const {
187  return track->GetCurrentStepNumber();
188  }
190  int pdgID() const {
191  G4ParticleDefinition* def = trackDef();
192  return def ? def->GetPDGEncoding() : 0;
193  }
195  const G4DynamicParticle* dynamic() const {
196  return track->GetDynamicParticle();
197  }
199  const G4PrimaryParticle* primary() const {
200  const G4DynamicParticle* d = track->GetDynamicParticle();
201  if ( d ) return d->GetPrimaryParticle();
202  return 0;
203  }
204  };
205 
206  } // End namespace Simulation
207 } // End namespace DD4hep
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.
Class of the Geant4 toolkit. See http://www-geant4.kek.jp/Reference.
Definition: Geant4Classes.h:60
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.
const G4Track * track
Reference to the track object.
const G4DynamicParticle * dynamic() const
Access the dynamic particle of the track object.
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.
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.