Geant4-11
Public Member Functions | Private Member Functions | Private Attributes
G4MaterialScanner Class Reference

#include <G4MaterialScanner.hh>

Public Member Functions

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

Private Member Functions

void DoScan ()
 
void RestoreUserActions ()
 
void StoreUserActions ()
 

Private Attributes

G4ThreeVector eyeDirection
 
G4ThreeVector eyePosition
 
G4int nPhi = 37
 
G4int nTheta = 91
 
G4double phiMin = 0.0
 
G4double phiSpan = 0.0
 
G4String regionName = "notDefined"
 
G4bool regionSensitive = false
 
G4EventManagertheEventManager = nullptr
 
G4UserEventActiontheMatScannerEventAction = nullptr
 
G4UserStackingActiontheMatScannerStackingAction = nullptr
 
G4MSSteppingActiontheMatScannerSteppingAction = nullptr
 
G4UserTrackingActiontheMatScannerTrackingAction = nullptr
 
G4MatScanMessengertheMessenger = nullptr
 
G4RayShootertheRayShooter = nullptr
 
G4RegiontheRegion = nullptr
 
G4double thetaMin = 0.0
 
G4double thetaSpan = 0.0
 
G4UserEventActiontheUserEventAction = nullptr
 
G4UserStackingActiontheUserStackingAction = nullptr
 
G4UserSteppingActiontheUserSteppingAction = nullptr
 
G4UserTrackingActiontheUserTrackingAction = nullptr
 

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{
55
56 eyePosition = G4ThreeVector(0., 0., 0.);
57 thetaSpan = 90. * deg;
58 phiSpan = 360. * deg;
59}
static constexpr double deg
Definition: G4SIunits.hh:132
CLHEP::Hep3Vector G4ThreeVector
static G4EventManager * GetEventManager()
G4MatScanMessenger * theMessenger
G4EventManager * theEventManager
G4ThreeVector eyePosition
G4RayShooter * theRayShooter

References deg, eyePosition, G4EventManager::GetEventManager(), phiSpan, theEventManager, theMessenger, theRayShooter, and thetaSpan.

◆ ~G4MaterialScanner()

G4MaterialScanner::~G4MaterialScanner ( )

Definition at line 62 of file G4MaterialScanner.cc.

63{
64 delete theRayShooter;
66 delete theMessenger;
67}
G4MSSteppingAction * theMatScannerSteppingAction

References theMatScannerSteppingAction, theMessenger, and theRayShooter.

Member Function Documentation

◆ DoScan()

void G4MaterialScanner::DoScan ( )
private

= G4Geantino::GeantinoDefinition()->GetProcessManager()->GetProcessList();

Definition at line 129 of file G4MaterialScanner.cc.

130{
131 // Confirm material table is updated
133
141
142 // Close geometry and set the application state
144 geomManager->OpenGeometry();
145 geomManager->CloseGeometry(1, 0);
146
147 G4ThreeVector center(0, 0, 0);
150 navigator->LocateGlobalPointAndSetup(center, 0, false);
151
153 theStateMan->SetNewState(G4State_GeomClosed);
154
155 // Event loop
156 G4int iEvent = 0;
157 for(G4int iTheta = 0; iTheta < nTheta; ++iTheta)
158 {
159 G4double theta = thetaMin;
160 if(iTheta > 0)
161 theta += G4double(iTheta) * thetaSpan / G4double(nTheta - 1);
162 G4double aveLength = 0.;
163 G4double aveX0 = 0.;
164 G4double aveLambda = 0.;
165 G4cout << G4endl;
166 G4cout
167 << " Theta(deg) Phi(deg) Length(mm) x0 lambda0"
168 << G4endl;
169 G4cout << G4endl;
170 for(G4int iPhi = 0; iPhi < nPhi; ++iPhi)
171 {
172 G4Event* anEvent = new G4Event(iEvent++);
173 G4double phi = phiMin;
174 if(iPhi > 0)
175 phi += G4double(iPhi) * phiSpan / G4double(nPhi - 1);
177 G4ThreeVector(std::cos(theta) * std::cos(phi),
178 std::cos(theta) * std::sin(phi), std::sin(theta));
185
186 G4cout << " " << std::setw(11) << theta / deg << " "
187 << std::setw(11) << phi / deg << " " << std::setw(11)
188 << length / mm << " " << std::setw(11) << x0 << " "
189 << std::setw(11) << lambda << G4endl;
190 aveLength += length / mm;
191 aveX0 += x0;
192 aveLambda += lambda;
193 }
194 if(nPhi > 1)
195 {
196 G4cout << G4endl;
197 G4cout << " ave. for theta = " << std::setw(11) << theta / deg << " : "
198 << std::setw(11) << aveLength / nPhi << " " << std::setw(11)
199 << aveX0 / nPhi << " " << std::setw(11) << aveLambda / nPhi
200 << G4endl;
201 }
202 }
203
204 theStateMan->SetNewState(G4State_Idle);
205 return;
206}
@ G4State_Idle
@ G4State_GeomClosed
static constexpr double mm
Definition: G4SIunits.hh:95
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void ProcessOneEvent(G4Event *anEvent)
static G4GeometryManager * GetInstance()
G4bool CloseGeometry(G4bool pOptimise=true, G4bool verbose=false, G4VPhysicalVolume *vol=nullptr)
void OpenGeometry(G4VPhysicalVolume *vol=nullptr)
G4double GetLambda0() const
G4double GetX0() const
G4double GetTotalStepLength() const
void Initialize(G4bool rSens, G4Region *reg)
G4ThreeVector eyeDirection
void Shoot(G4Event *evt, G4ThreeVector vtx, G4ThreeVector direc)
Definition: G4RayShooter.cc:57
static G4RunManagerKernel * GetRunManagerKernel()
static G4StateManager * GetStateManager()
G4bool SetNewState(const G4ApplicationState &requestedState)
static G4TransportationManager * GetTransportationManager()
G4Navigator * GetNavigatorForTracking() const

References G4GeometryManager::CloseGeometry(), deg, eyeDirection, eyePosition, G4cout, G4endl, G4State_GeomClosed, G4State_Idle, G4GeometryManager::GetInstance(), G4MSSteppingAction::GetLambda0(), G4TransportationManager::GetNavigatorForTracking(), G4RunManagerKernel::GetRunManagerKernel(), G4StateManager::GetStateManager(), G4MSSteppingAction::GetTotalStepLength(), G4TransportationManager::GetTransportationManager(), G4MSSteppingAction::GetX0(), G4MSSteppingAction::Initialize(), G4InuclParticleNames::lambda, mm, write_gdml::navigator, nPhi, nTheta, G4GeometryManager::OpenGeometry(), phiMin, phiSpan, G4EventManager::ProcessOneEvent(), regionSensitive, G4StateManager::SetNewState(), G4RayShooter::Shoot(), theEventManager, theMatScannerSteppingAction, theRayShooter, theRegion, thetaMin, thetaSpan, and G4RunManagerKernel::UpdateRegion().

Referenced by Scan().

◆ GetEyePosition()

G4ThreeVector G4MaterialScanner::GetEyePosition ( ) const
inline

Definition at line 64 of file G4MaterialScanner.hh.

64{ return eyePosition; }

References eyePosition.

Referenced by G4MatScanMessenger::GetCurrentValue().

◆ GetNPhi()

G4int G4MaterialScanner::GetNPhi ( ) const
inline

Definition at line 72 of file G4MaterialScanner.hh.

72{ return nPhi; }

References 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; }

References 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; }

References 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; }

References 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; }

References regionName.

Referenced by G4MatScanMessenger::GetCurrentValue().

◆ GetRegionSensitive()

G4bool G4MaterialScanner::GetRegionSensitive ( ) const
inline

Definition at line 78 of file G4MaterialScanner.hh.

78{ return regionSensitive; }

References regionSensitive.

Referenced by G4MatScanMessenger::GetCurrentValue().

◆ GetThetaMin()

G4double G4MaterialScanner::GetThetaMin ( ) const
inline

Definition at line 68 of file G4MaterialScanner.hh.

68{ return thetaMin; }

References 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; }

References thetaSpan.

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

◆ RestoreUserActions()

void G4MaterialScanner::RestoreUserActions ( )
private

Definition at line 114 of file G4MaterialScanner.cc.

115{
120
122 if(theSDMan != nullptr)
123 {
124 theSDMan->Activate("/", true);
125 }
126}
void SetUserAction(G4UserEventAction *userAction)
G4UserTrackingAction * theUserTrackingAction
G4UserStackingAction * theUserStackingAction
G4UserSteppingAction * theUserSteppingAction
G4UserEventAction * theUserEventAction
void Activate(G4String dName, G4bool activeFlag)
Definition: G4SDManager.cc:125
static G4SDManager * GetSDMpointerIfExist()
Definition: G4SDManager.cc:47

References G4SDManager::Activate(), G4SDManager::GetSDMpointerIfExist(), G4EventManager::SetUserAction(), theEventManager, theUserEventAction, theUserStackingAction, theUserSteppingAction, and theUserTrackingAction.

Referenced by Scan().

◆ 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 {
83 }
85 DoScan();
87}
G4ApplicationState
G4GLOB_DLL std::ostream G4cerr
const G4ApplicationState & GetCurrentState() const

References DoScan(), G4cerr, G4endl, G4State_Idle, G4StateManager::GetCurrentState(), G4StateManager::GetStateManager(), RestoreUserActions(), StoreUserActions(), and theMatScannerSteppingAction.

Referenced by G4MatScanMessenger::SetNewValue().

◆ SetEyePosition()

void G4MaterialScanner::SetEyePosition ( const G4ThreeVector val)
inline

Definition at line 63 of file G4MaterialScanner.hh.

63{ eyePosition = val; }

References eyePosition.

Referenced by G4MatScanMessenger::SetNewValue().

◆ SetNPhi()

void G4MaterialScanner::SetNPhi ( G4int  val)
inline

Definition at line 71 of file G4MaterialScanner.hh.

71{ nPhi = val; }

References nPhi.

Referenced by G4MatScanMessenger::SetNewValue().

◆ SetNTheta()

void G4MaterialScanner::SetNTheta ( G4int  val)
inline

Definition at line 65 of file G4MaterialScanner.hh.

65{ nTheta = val; }

References nTheta.

Referenced by G4MatScanMessenger::SetNewValue().

◆ SetPhiMin()

void G4MaterialScanner::SetPhiMin ( G4double  val)
inline

Definition at line 73 of file G4MaterialScanner.hh.

73{ phiMin = val; }

References phiMin.

Referenced by G4MatScanMessenger::SetNewValue().

◆ SetPhiSpan()

void G4MaterialScanner::SetPhiSpan ( G4double  val)
inline

Definition at line 75 of file G4MaterialScanner.hh.

75{ phiSpan = val; }

References phiSpan.

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

References G4cerr, G4endl, G4RegionStore::GetInstance(), G4RegionStore::GetRegion(), regionName, and theRegion.

Referenced by G4MatScanMessenger::SetNewValue().

◆ SetRegionSensitive()

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

Definition at line 77 of file G4MaterialScanner.hh.

77{ regionSensitive = val; }

References regionSensitive.

Referenced by G4MatScanMessenger::SetNewValue().

◆ SetThetaMin()

void G4MaterialScanner::SetThetaMin ( G4double  val)
inline

Definition at line 67 of file G4MaterialScanner.hh.

67{ thetaMin = val; }

References thetaMin.

Referenced by G4MatScanMessenger::SetNewValue().

◆ SetThetaSpan()

void G4MaterialScanner::SetThetaSpan ( G4double  val)
inline

Definition at line 69 of file G4MaterialScanner.hh.

69{ thetaSpan = val; }

References thetaSpan.

Referenced by G4MatScanMessenger::SetNewValue().

◆ StoreUserActions()

void G4MaterialScanner::StoreUserActions ( )
private

Definition at line 90 of file G4MaterialScanner.cc.

91{
96
101
103 if(theSDMan != nullptr)
104 {
105 theSDMan->Activate("/", false);
106 }
107
109 theGeomMan->OpenGeometry();
110 theGeomMan->CloseGeometry(true);
111}
G4UserTrackingAction * GetUserTrackingAction()
G4UserSteppingAction * GetUserSteppingAction()
G4UserStackingAction * GetUserStackingAction()
G4UserEventAction * GetUserEventAction()
G4UserStackingAction * theMatScannerStackingAction
G4UserTrackingAction * theMatScannerTrackingAction
G4UserEventAction * theMatScannerEventAction

References G4SDManager::Activate(), G4GeometryManager::CloseGeometry(), G4GeometryManager::GetInstance(), G4SDManager::GetSDMpointerIfExist(), G4EventManager::GetUserEventAction(), G4EventManager::GetUserStackingAction(), G4EventManager::GetUserSteppingAction(), G4EventManager::GetUserTrackingAction(), G4GeometryManager::OpenGeometry(), G4EventManager::SetUserAction(), theEventManager, theMatScannerEventAction, theMatScannerStackingAction, theMatScannerSteppingAction, theMatScannerTrackingAction, theUserEventAction, theUserStackingAction, theUserSteppingAction, and theUserTrackingAction.

Referenced by Scan().

Field Documentation

◆ eyeDirection

G4ThreeVector G4MaterialScanner::eyeDirection
private

Definition at line 116 of file G4MaterialScanner.hh.

Referenced by DoScan().

◆ eyePosition

G4ThreeVector G4MaterialScanner::eyePosition
private

Definition at line 108 of file G4MaterialScanner.hh.

Referenced by DoScan(), G4MaterialScanner(), GetEyePosition(), and SetEyePosition().

◆ nPhi

G4int G4MaterialScanner::nPhi = 37
private

Definition at line 112 of file G4MaterialScanner.hh.

Referenced by DoScan(), GetNPhi(), and SetNPhi().

◆ nTheta

G4int G4MaterialScanner::nTheta = 91
private

Definition at line 109 of file G4MaterialScanner.hh.

Referenced by DoScan(), GetNTheta(), and SetNTheta().

◆ phiMin

G4double G4MaterialScanner::phiMin = 0.0
private

Definition at line 113 of file G4MaterialScanner.hh.

Referenced by DoScan(), GetPhiMin(), and SetPhiMin().

◆ phiSpan

G4double G4MaterialScanner::phiSpan = 0.0
private

Definition at line 114 of file G4MaterialScanner.hh.

Referenced by DoScan(), G4MaterialScanner(), GetPhiSpan(), and SetPhiSpan().

◆ regionName

G4String G4MaterialScanner::regionName = "notDefined"
private

Definition at line 119 of file G4MaterialScanner.hh.

Referenced by GetRegionName(), and SetRegionName().

◆ regionSensitive

G4bool G4MaterialScanner::regionSensitive = false
private

Definition at line 118 of file G4MaterialScanner.hh.

Referenced by DoScan(), GetRegionSensitive(), and SetRegionSensitive().

◆ theEventManager

G4EventManager* G4MaterialScanner::theEventManager = nullptr
private

◆ theMatScannerEventAction

G4UserEventAction* G4MaterialScanner::theMatScannerEventAction = nullptr
private

Definition at line 103 of file G4MaterialScanner.hh.

Referenced by StoreUserActions().

◆ theMatScannerStackingAction

G4UserStackingAction* G4MaterialScanner::theMatScannerStackingAction = nullptr
private

Definition at line 104 of file G4MaterialScanner.hh.

Referenced by StoreUserActions().

◆ theMatScannerSteppingAction

G4MSSteppingAction* G4MaterialScanner::theMatScannerSteppingAction = nullptr
private

Definition at line 106 of file G4MaterialScanner.hh.

Referenced by DoScan(), Scan(), StoreUserActions(), and ~G4MaterialScanner().

◆ theMatScannerTrackingAction

G4UserTrackingAction* G4MaterialScanner::theMatScannerTrackingAction = nullptr
private

Definition at line 105 of file G4MaterialScanner.hh.

Referenced by StoreUserActions().

◆ theMessenger

G4MatScanMessenger* G4MaterialScanner::theMessenger = nullptr
private

Definition at line 94 of file G4MaterialScanner.hh.

Referenced by G4MaterialScanner(), and ~G4MaterialScanner().

◆ theRayShooter

G4RayShooter* G4MaterialScanner::theRayShooter = nullptr
private

Definition at line 93 of file G4MaterialScanner.hh.

Referenced by DoScan(), G4MaterialScanner(), and ~G4MaterialScanner().

◆ theRegion

G4Region* G4MaterialScanner::theRegion = nullptr
private

Definition at line 120 of file G4MaterialScanner.hh.

Referenced by DoScan(), and SetRegionName().

◆ thetaMin

G4double G4MaterialScanner::thetaMin = 0.0
private

Definition at line 110 of file G4MaterialScanner.hh.

Referenced by DoScan(), GetThetaMin(), and SetThetaMin().

◆ thetaSpan

G4double G4MaterialScanner::thetaSpan = 0.0
private

Definition at line 111 of file G4MaterialScanner.hh.

Referenced by DoScan(), G4MaterialScanner(), GetThetaSpan(), and SetThetaSpan().

◆ theUserEventAction

G4UserEventAction* G4MaterialScanner::theUserEventAction = nullptr
private

Definition at line 98 of file G4MaterialScanner.hh.

Referenced by RestoreUserActions(), and StoreUserActions().

◆ theUserStackingAction

G4UserStackingAction* G4MaterialScanner::theUserStackingAction = nullptr
private

Definition at line 99 of file G4MaterialScanner.hh.

Referenced by RestoreUserActions(), and StoreUserActions().

◆ theUserSteppingAction

G4UserSteppingAction* G4MaterialScanner::theUserSteppingAction = nullptr
private

Definition at line 101 of file G4MaterialScanner.hh.

Referenced by RestoreUserActions(), and StoreUserActions().

◆ theUserTrackingAction

G4UserTrackingAction* G4MaterialScanner::theUserTrackingAction = nullptr
private

Definition at line 100 of file G4MaterialScanner.hh.

Referenced by RestoreUserActions(), and StoreUserActions().


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