|
DD4hep
1.30.0
Detector Description Toolkit for High Energy Physics
|
Go to the documentation of this file.
13 #ifndef DDG4_GEANT4STEPHANDLER_H
14 #define DDG4_GEANT4STEPHANDLER_H
21 #include <G4StepPoint.hh>
22 #include <G4VTouchable.hh>
23 #include <G4VSensitiveDetector.hh>
24 #include <G4EmSaturation.hh>
25 #include <G4Version.hh>
34 class Geant4StepHandler;
56 :
Geant4HitHandler(s->GetTrack(), (s->GetPreStepPoint()->GetTouchableHandle())()),
57 step(s),
pre(s->GetPreStepPoint()),
post(s->GetPostStepPoint())
71 static const char*
stepStatus(G4StepStatus status);
81 return step->GetTotalEnergyDeposit();
85 const G4ThreeVector& p =
pre->GetPosition();
86 return Position(p.x(), p.y(), p.z());
90 return pre->GetPosition();
94 const G4ThreeVector& p =
post->GetPosition();
95 return Position(p.x(), p.y(), p.z());
99 return post->GetPosition();
103 const G4ThreeVector& p1 =
pre->GetPosition();
104 const G4ThreeVector& p2 =
post->GetPosition();
105 G4ThreeVector r = 0.5*(p2+p1);
110 const G4ThreeVector& p1 =
pre->GetPosition();
111 const G4ThreeVector& p2 =
post->GetPosition();
112 G4ThreeVector r = 0.5*(p2+p1);
117 const G4ThreeVector& p1 =
pre->GetPosition();
118 const G4ThreeVector& p2 =
post->GetPosition();
119 G4ThreeVector r = (p2-p1);
128 const G4ThreeVector& p =
pre->GetMomentum();
129 return Momentum(p.x(), p.y(), p.z());
132 const G4ThreeVector& p =
post->GetMomentum();
133 return Momentum(p.x(), p.y(), p.z());
136 return step->GetTotalEnergyDeposit();
139 return step->GetStepLength();
142 return pre->GetTouchable();
145 return post->GetTouchable();
147 const char*
volName(
const G4StepPoint* p,
const char* undefined =
"")
const {
148 G4VPhysicalVolume*
v =
volume(p);
149 return v ?
v->GetName().c_str() : undefined;
151 G4VPhysicalVolume*
volume(
const G4StepPoint* p)
const {
152 return p->GetTouchableHandle()->GetVolume();
154 G4VSolid*
solid(
const G4StepPoint* p)
const {
155 return p->GetTouchableHandle()->GetSolid();
157 G4VPhysicalVolume*
physvol(
const G4StepPoint* p)
const {
158 return p->GetPhysicalVolume();
160 G4LogicalVolume*
logvol(
const G4StepPoint* p)
const {
161 G4VPhysicalVolume* pv =
physvol(p);
162 return pv ? pv->GetLogicalVolume() : 0;
165 G4LogicalVolume* lv =
logvol(p);
166 return lv ? lv->GetSensitiveDetector() : 0;
168 std::string
sdName(
const G4StepPoint* p,
const std::string& undefined =
"")
const {
170 return s ? s->GetName().c_str() : undefined;
188 return step->IsFirstStepInVolume();
191 return step->IsLastStepInVolume();
204 #endif // DDG4_GEANT4STEPHANDLER_H
Geant4StepHandler(const Geant4StepHandler ©)=delete
No copy constructor.
Position prePos() const
Returns the pre-step position.
Geant4StepHandler & operator=(Geant4StepHandler &©)=delete
Move operator inhibited. Should not be copied.
const char * volName(const G4StepPoint *p, const char *undefined="") const
Geant4StepHandler()=delete
Inhibit default constructor.
bool isSensitive(const G4StepPoint *point) const
const G4VTouchable * postTouchable() const
const G4ThreeVector & prePosG4() const
Returns the pre-step position as a G4ThreeVector.
G4VPhysicalVolume * volume(const G4StepPoint *p) const
G4VPhysicalVolume * preVolume() const
bool lastInVolume() const
const char * postStepStatus() const
Returns the post-step status in form of a string.
Helper class to ease the extraction of information from a G4Step object.
G4ThreeVector avgPositionG4() const
Average position vector of the step.
Geant4StepHandler(Geant4StepHandler &©)=delete
No move constructor.
double stepLength() const
Helper class to ease the extraction of information from a G4Hit object.
Class of the Geant4 toolkit. See http://www-geant4.kek.jp/Reference.
G4VPhysicalVolume * physvol(const G4StepPoint *p) const
G4VSolid * solid(const G4StepPoint *p) const
Geant4StepHandler & operator=(const Geant4StepHandler ©)=delete
Assignment operator inhibited. Should not be copied.
bool isSensitive(const G4LogicalVolume *lv) const
G4VPhysicalVolume * postVolume() const
const G4ThreeVector & postPosG4() const
Returns the post-step position as a G4ThreeVector.
Position postPos() const
Returns the post-step position.
double birkAttenuation() const
Apply BirksLaw.
Position direction() const
Direction calculated from the post- and pre-position ofthe step.
ROOT::Math::XYZVector Position
G4VSensitiveDetector * preSD() const
static const char * stepStatus(G4StepStatus status)
Returns the step status (argument) in form of a string.
Namespace for the AIDA detector description toolkit.
std::string sdName(const G4StepPoint *p, const std::string &undefined="") const
void doApplyBirksLaw(void)
Set applyBirksLaw to ture.
G4VSensitiveDetector * postSD() const
const char * preStepStatus() const
Returns the pre-step status in form of a string.
Geant4StepHandler(const G4Step *s)
Initializing constructor.
Position avgPosition() const
Average position vector of the step.
const G4VTouchable * preTouchable() const
G4VSensitiveDetector * sd(const G4StepPoint *p) const
double totalEnergy() const
Returns total energy deposit.
bool firstInVolume() const
G4LogicalVolume * logvol(const G4StepPoint *p) const