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

#include <G4RTSteppingAction.hh>

+ Inheritance diagram for G4RTSteppingAction:

Public Member Functions

 G4RTSteppingAction ()
 
virtual ~G4RTSteppingAction ()
 
virtual void UserSteppingAction (const G4Step *)
 
- Public Member Functions inherited from G4UserSteppingAction
 G4UserSteppingAction ()
 
virtual ~G4UserSteppingAction ()
 
virtual void SetSteppingManagerPointer (G4SteppingManager *pValue)
 
virtual void UserSteppingAction (const G4Step *)
 

Static Public Member Functions

static void SetIgnoreTransparency (G4bool val)
 
static G4bool GetIgnoreTransparency ()
 

Additional Inherited Members

- Protected Attributes inherited from G4UserSteppingAction
G4SteppingManagerfpSteppingManager = nullptr
 

Detailed Description

Definition at line 51 of file G4RTSteppingAction.hh.

Constructor & Destructor Documentation

◆ G4RTSteppingAction()

G4RTSteppingAction::G4RTSteppingAction ( )

Definition at line 45 of file G4RTSteppingAction.cc.

46{}

◆ ~G4RTSteppingAction()

virtual G4RTSteppingAction::~G4RTSteppingAction ( )
inlinevirtual

Definition at line 55 of file G4RTSteppingAction.hh.

55{;}

Member Function Documentation

◆ GetIgnoreTransparency()

G4bool G4RTSteppingAction::GetIgnoreTransparency ( )
static

Definition at line 99 of file G4RTSteppingAction.cc.

100{ return ignoreTransparency; }

Referenced by G4RTMessenger::GetCurrentValue().

◆ SetIgnoreTransparency()

void G4RTSteppingAction::SetIgnoreTransparency ( G4bool  val)
static

Definition at line 97 of file G4RTSteppingAction.cc.

98{ ignoreTransparency = val; }

Referenced by G4RTMessenger::SetNewValue().

◆ UserSteppingAction()

void G4RTSteppingAction::UserSteppingAction ( const G4Step aStep)
virtual

Reimplemented from G4UserSteppingAction.

Definition at line 48 of file G4RTSteppingAction.cc.

49{
50 // Stop the ray particle if it reaches to the coloured volume with its alpha = 1
51 // or coloured volume and the user request to ignore transparency
52
53 G4StepPoint* postStepPoint = aStep -> GetPostStepPoint();
54
55 G4VPhysicalVolume* postVolume_phys=postStepPoint->GetPhysicalVolume();
56 if(!postVolume_phys) return; // Reach to out of the world
57
59 G4RayTracerSceneHandler* sceneHandler =
60 static_cast<G4RayTracerSceneHandler*>(visManager->GetCurrentSceneHandler());
61
62 // Make a path from the touchable
63 const G4VTouchable* postTouchable = postStepPoint->GetTouchable();
64 G4int postDepth = postTouchable->GetHistoryDepth();
65 G4ModelingParameters::PVPointerCopyNoPath localPostPVPointerCopyNoPath;
66 for (G4int i = postDepth; i >= 0; --i) {
67 localPostPVPointerCopyNoPath.push_back
69 (postTouchable->GetVolume(i),postTouchable->GetCopyNumber(i)));
70 }
71
72 // Pick up the vis atts, if any, from the scene handler
73 const auto& sceneVisAttsMap = sceneHandler->GetSceneVisAttsMap();
74 auto postIterator = sceneVisAttsMap.find(localPostPVPointerCopyNoPath);
75 const G4VisAttributes* postVisAtts;
76 if (postIterator != sceneVisAttsMap.end()) {
77 postVisAtts = &postIterator->second;
78 } else {
79 postVisAtts = 0;
80 }
81
82 if((!postVisAtts) || (!(postVisAtts->IsVisible()))) return; // Invisible volume, continue tracking
83
84 if((postVisAtts->IsForceDrawingStyle())
85 &&(postVisAtts->GetForcedDrawingStyle()==G4VisAttributes::wireframe)) return;
86 // Wire frame volume, continue tracking
87
88 G4double postAlpha=(postVisAtts->GetColour()).GetAlpha();
89
90 if(postAlpha==1.0 || ignoreTransparency) // Stop stepping
91 {
92 G4Track* currentTrack = aStep -> GetTrack();
93 currentTrack -> SetTrackStatus(fStopAndKill);
94 }
95}
@ fStopAndKill
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
std::vector< PVPointerCopyNo > PVPointerCopyNoPath
const std::map< G4ModelingParameters::PVPointerCopyNoPath, G4VisAttributes, PathLessThan > & GetSceneVisAttsMap() const
const G4VTouchable * GetTouchable() const
G4VPhysicalVolume * GetPhysicalVolume() const
G4int GetCopyNumber(G4int depth=0) const
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:41
virtual G4int GetHistoryDepth() const
Definition: G4VTouchable.cc:76
const G4Colour & GetColour() const
G4bool IsVisible() const
ForcedDrawingStyle GetForcedDrawingStyle() const
G4bool IsForceDrawingStyle() const
G4VSceneHandler * GetCurrentSceneHandler() const
static G4VisManager * GetInstance()

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