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

#include <G4MaterialScanner.hh>

Public Member Functions

 G4MaterialScanner ()
 
 ~G4MaterialScanner ()
 
void Scan ()
 
void SetEyePosition (const G4ThreeVector &val)
 
G4ThreeVector GetEyePosition () const
 
void SetNTheta (G4int val)
 
G4int GetNTheta () const
 
void SetThetaMin (G4double val)
 
G4double GetThetaMin () const
 
void SetThetaSpan (G4double val)
 
G4double GetThetaSpan () const
 
void SetNPhi (G4int val)
 
G4int GetNPhi () const
 
void SetPhiMin (G4double val)
 
G4double GetPhiMin () const
 
void SetPhiSpan (G4double val)
 
G4double GetPhiSpan () const
 
void SetRegionSensitive (G4bool val=true)
 
G4bool GetRegionSensitive () const
 
G4bool SetRegionName (const G4String &val)
 
const G4StringGetRegionName () const
 

Detailed Description

Definition at line 52 of file G4MaterialScanner.hh.

Constructor & Destructor Documentation

◆ G4MaterialScanner()

G4MaterialScanner::G4MaterialScanner ( )

Definition at line 50 of file G4MaterialScanner.cc.

51{
52 theRayShooter = new G4RayShooter();
53 theMessenger = new G4MatScanMessenger(this);
54 theEventManager = G4EventManager::GetEventManager();
55
56 eyePosition = G4ThreeVector(0., 0., 0.);
57 thetaSpan = 90. * deg;
58 phiSpan = 360. * deg;
59}
CLHEP::Hep3Vector G4ThreeVector
static G4EventManager * GetEventManager()

◆ ~G4MaterialScanner()

G4MaterialScanner::~G4MaterialScanner ( )

Definition at line 62 of file G4MaterialScanner.cc.

63{
64 delete theRayShooter;
65 delete theMatScannerSteppingAction;
66 delete theMessenger;
67}

Member Function Documentation

◆ GetEyePosition()

G4ThreeVector G4MaterialScanner::GetEyePosition ( ) const
inline

Definition at line 64 of file G4MaterialScanner.hh.

64{ return eyePosition; }

Referenced by G4MatScanMessenger::GetCurrentValue().

◆ GetNPhi()

G4int G4MaterialScanner::GetNPhi ( ) const
inline

Definition at line 72 of file G4MaterialScanner.hh.

72{ return nPhi; }

Referenced by G4MatScanMessenger::GetCurrentValue(), and G4MatScanMessenger::SetNewValue().

◆ GetNTheta()

G4int G4MaterialScanner::GetNTheta ( ) const
inline

Definition at line 66 of file G4MaterialScanner.hh.

66{ return nTheta; }

Referenced by G4MatScanMessenger::GetCurrentValue(), and G4MatScanMessenger::SetNewValue().

◆ GetPhiMin()

G4double G4MaterialScanner::GetPhiMin ( ) const
inline

Definition at line 74 of file G4MaterialScanner.hh.

74{ return phiMin; }

Referenced by G4MatScanMessenger::GetCurrentValue(), and G4MatScanMessenger::SetNewValue().

◆ GetPhiSpan()

G4double G4MaterialScanner::GetPhiSpan ( ) const
inline

Definition at line 76 of file G4MaterialScanner.hh.

76{ return phiSpan; }

Referenced by G4MatScanMessenger::GetCurrentValue(), and G4MatScanMessenger::SetNewValue().

◆ GetRegionName()

const G4String & G4MaterialScanner::GetRegionName ( ) const
inline

Definition at line 80 of file G4MaterialScanner.hh.

80{ return regionName; }

Referenced by G4MatScanMessenger::GetCurrentValue().

◆ GetRegionSensitive()

G4bool G4MaterialScanner::GetRegionSensitive ( ) const
inline

Definition at line 78 of file G4MaterialScanner.hh.

78{ return regionSensitive; }

Referenced by G4MatScanMessenger::GetCurrentValue().

◆ GetThetaMin()

G4double G4MaterialScanner::GetThetaMin ( ) const
inline

Definition at line 68 of file G4MaterialScanner.hh.

68{ return thetaMin; }

Referenced by G4MatScanMessenger::GetCurrentValue(), and G4MatScanMessenger::SetNewValue().

◆ GetThetaSpan()

G4double G4MaterialScanner::GetThetaSpan ( ) const
inline

Definition at line 70 of file G4MaterialScanner.hh.

70{ return thetaSpan; }

Referenced by G4MatScanMessenger::GetCurrentValue(), and G4MatScanMessenger::SetNewValue().

◆ Scan()

void G4MaterialScanner::Scan ( )

Definition at line 70 of file G4MaterialScanner.cc.

71{
73 G4ApplicationState currentState = theStateMan->GetCurrentState();
74 if(currentState != G4State_Idle)
75 {
76 G4cerr << "Illegal application state - Scan() ignored." << G4endl;
77 return;
78 }
79
80 if(theMatScannerSteppingAction == nullptr)
81 {
82 theMatScannerSteppingAction = new G4MSSteppingAction();
83 }
84 StoreUserActions();
85 DoScan();
86 RestoreUserActions();
87}
G4ApplicationState
@ G4State_Idle
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition: G4ios.hh:57
const G4ApplicationState & GetCurrentState() const
static G4StateManager * GetStateManager()

Referenced by G4MatScanMessenger::SetNewValue().

◆ SetEyePosition()

void G4MaterialScanner::SetEyePosition ( const G4ThreeVector val)
inline

Definition at line 63 of file G4MaterialScanner.hh.

63{ eyePosition = val; }

Referenced by G4MatScanMessenger::SetNewValue().

◆ SetNPhi()

void G4MaterialScanner::SetNPhi ( G4int  val)
inline

Definition at line 71 of file G4MaterialScanner.hh.

71{ nPhi = val; }

Referenced by G4MatScanMessenger::SetNewValue().

◆ SetNTheta()

void G4MaterialScanner::SetNTheta ( G4int  val)
inline

Definition at line 65 of file G4MaterialScanner.hh.

65{ nTheta = val; }

Referenced by G4MatScanMessenger::SetNewValue().

◆ SetPhiMin()

void G4MaterialScanner::SetPhiMin ( G4double  val)
inline

Definition at line 73 of file G4MaterialScanner.hh.

73{ phiMin = val; }

Referenced by G4MatScanMessenger::SetNewValue().

◆ SetPhiSpan()

void G4MaterialScanner::SetPhiSpan ( G4double  val)
inline

Definition at line 75 of file G4MaterialScanner.hh.

75{ phiSpan = val; }

Referenced by G4MatScanMessenger::SetNewValue().

◆ SetRegionName()

G4bool G4MaterialScanner::SetRegionName ( const G4String val)

Definition at line 209 of file G4MaterialScanner.cc.

210{
212 if(aRegion != nullptr)
213 {
214 theRegion = aRegion;
215 regionName = val;
216 return true;
217 }
218 else
219 {
220 G4cerr << "Region <" << val << "> not found. Command ignored." << G4endl;
221 G4cerr << "Defined regions are : " << G4endl;
222 for(std::size_t i = 0; i < G4RegionStore::GetInstance()->size(); ++i)
223 {
224 G4cerr << " " << (*(G4RegionStore::GetInstance()))[i]->GetName();
225 }
226 G4cerr << G4endl;
227 return false;
228 }
229}
static G4RegionStore * GetInstance()
G4Region * GetRegion(const G4String &name, G4bool verbose=true) const

Referenced by G4MatScanMessenger::SetNewValue().

◆ SetRegionSensitive()

void G4MaterialScanner::SetRegionSensitive ( G4bool  val = true)
inline

Definition at line 77 of file G4MaterialScanner.hh.

77{ regionSensitive = val; }

Referenced by G4MatScanMessenger::SetNewValue().

◆ SetThetaMin()

void G4MaterialScanner::SetThetaMin ( G4double  val)
inline

Definition at line 67 of file G4MaterialScanner.hh.

67{ thetaMin = val; }

Referenced by G4MatScanMessenger::SetNewValue().

◆ SetThetaSpan()

void G4MaterialScanner::SetThetaSpan ( G4double  val)
inline

Definition at line 69 of file G4MaterialScanner.hh.

69{ thetaSpan = val; }

Referenced by G4MatScanMessenger::SetNewValue().


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