8#include "G4LogicalVolume.hh"
9#include "G4VPhysicalVolume.hh"
10#include "G4TouchableHistory.hh"
11#include "G4VProcess.hh"
16 assert (theStep!=
nullptr);
17 return theStep->GetPreStepPoint()->GetPosition();
21 assert (theStep!=
nullptr);
22 return theStep->GetPostStepPoint()->GetPosition();
26 assert (theStep!=
nullptr);
27 return theStep->GetTrack()->GetDefinition()->GetParticleName();
31 assert (theStep!=
nullptr);
32 return theStep->GetTrack()->GetDefinition()->GetPDGEncoding();
36 assert (theStep!=
nullptr);
37 return theStep->GetTotalEnergyDeposit();
49 const G4TouchableHistory *history
50 =
static_cast<const G4TouchableHistory*
>(theStep->GetPreStepPoint()->GetTouchable());
53 return history->GetVolume(std::abs(iLevel));
55 const int nLev=history->GetHistoryDepth();
56 return history->GetVolume(nLev-iLevel);
68 const G4TouchableHistory *history
69 =
static_cast<const G4TouchableHistory*
>(theStep->GetPostStepPoint()->GetTouchable());
72 return history->GetVolume(std::abs(iLevel));
74 const int nLev=history->GetHistoryDepth();
75 return history->GetVolume(nLev-iLevel);
79 return static_cast<const G4TouchableHistory*
>(theStep->GetPreStepPoint()->
80 GetTouchable())->GetHistoryDepth();
84 return static_cast<const G4TouchableHistory*
>(theStep->GetPostStepPoint()->
85 GetTouchable())->GetHistoryDepth();
89 return theStep->GetPostStepPoint()->GetProcessDefinedStep();
99 return getProcess(theStep)->GetProcessSubType();
helper functions to avoid having to play with the G4Step to retrieve relevant quantities.
std::string getPreStepLogicalVolumeName(const G4Step *theStep, int iLevel=0)
TODO.
int particlePDGCode(const G4Step *theStep)
TODO.
double depositedEnergy(const G4Step *theStep)
TODO.
int postStepBranchDepth(const G4Step *theStep)
TODO.
G4VPhysicalVolume * getPreStepPhysicalVolume(const G4Step *theStep, int iLevel=0)
TODO.
G4ThreeVector preStepPosition(const G4Step *theStep)
TODO.
int preStepBranchDepth(const G4Step *theStep)
TODO.
const G4VProcess * getProcess(const G4Step *theStep)
TODO.
std::string particleName(const G4Step *theStep)
TODO.
G4ThreeVector postStepPosition(const G4Step *theStep)
TODO.
std::string getPostStepLogicalVolumeName(const G4Step *theStep, int iLevel=0)
TODO.
G4int getProcessSubType(const G4Step *theStep)
TODO.
G4LogicalVolume * getPreStepLogicalVolume(const G4Step *theStep, int iLevel=0)
TODO.
G4LogicalVolume * getPostStepLogicalVolume(const G4Step *theStep, int iLevel=0)
TODO.
std::string getProcessName(const G4Step *theStep)
TODO.
G4VPhysicalVolume * getPostStepPhysicalVolume(const G4Step *theStep, int iLevel=0)
TODO.