Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4MSSteppingAction Class Reference

#include <G4MSSteppingAction.hh>

+ Inheritance diagram for G4MSSteppingAction:

Public Member Functions

 G4MSSteppingAction ()
 
virtual ~G4MSSteppingAction ()
 
void Initialize (G4bool rSens, G4Region *reg)
 
virtual void UserSteppingAction (const G4Step *)
 
G4double GetTotalStepLength () const
 
G4double GetX0 () const
 
G4double GetLambda0 () const
 
- Public Member Functions inherited from G4UserSteppingAction
 G4UserSteppingAction ()
 
virtual ~G4UserSteppingAction ()
 
void SetSteppingManagerPointer (G4SteppingManager *pValue)
 
virtual void UserSteppingAction (const G4Step *)
 

Additional Inherited Members

- Protected Attributes inherited from G4UserSteppingAction
G4SteppingManagerfpSteppingManager
 

Detailed Description

Definition at line 47 of file G4MSSteppingAction.hh.

Constructor & Destructor Documentation

◆ G4MSSteppingAction()

G4MSSteppingAction::G4MSSteppingAction ( )

Definition at line 41 of file G4MSSteppingAction.cc.

42{
43 Initialize(false,0);
44}
void Initialize(G4bool rSens, G4Region *reg)

◆ ~G4MSSteppingAction()

G4MSSteppingAction::~G4MSSteppingAction ( )
virtual

Definition at line 46 of file G4MSSteppingAction.cc.

47{;}

Member Function Documentation

◆ GetLambda0()

G4double G4MSSteppingAction::GetLambda0 ( ) const
inline

Definition at line 66 of file G4MSSteppingAction.hh.

66{ return lambda; }

◆ GetTotalStepLength()

G4double G4MSSteppingAction::GetTotalStepLength ( ) const
inline

Definition at line 64 of file G4MSSteppingAction.hh.

64{ return length; }

◆ GetX0()

G4double G4MSSteppingAction::GetX0 ( ) const
inline

Definition at line 65 of file G4MSSteppingAction.hh.

65{ return x0; }

◆ Initialize()

void G4MSSteppingAction::Initialize ( G4bool  rSens,
G4Region reg 
)

Definition at line 49 of file G4MSSteppingAction.cc.

50{
51 regionSensitive = rSens;
52 theRegion = reg;
53 length = 0.;
54 x0 = 0.;
55 lambda = 0.;
56}

Referenced by G4MSSteppingAction().

◆ UserSteppingAction()

void G4MSSteppingAction::UserSteppingAction ( const G4Step aStep)
virtual

Reimplemented from G4UserSteppingAction.

Definition at line 58 of file G4MSSteppingAction.cc.

59{
60 G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
61 G4Region* region = preStepPoint->GetPhysicalVolume()->GetLogicalVolume()->GetRegion();
62
63 if(regionSensitive && (region!=theRegion)) return;
64
65 G4double stlen = aStep->GetStepLength();
66 G4Material* material = preStepPoint->GetMaterial();
67 length += stlen;
68 x0 += stlen/(material->GetRadlen());
69 lambda += stlen/(material->GetNuclearInterLength());
70}
double G4double
Definition: G4Types.hh:64
G4Region * GetRegion() const
G4double GetRadlen() const
Definition: G4Material.hh:219
G4double GetNuclearInterLength() const
Definition: G4Material.hh:222
G4Material * GetMaterial() const
G4VPhysicalVolume * GetPhysicalVolume() const
G4StepPoint * GetPreStepPoint() const
G4double GetStepLength() const
G4LogicalVolume * GetLogicalVolume() const

The documentation for this class was generated from the following files: