Geant4  10.02.p03
G4MaterialScanner Class Reference

#include <G4MaterialScanner.hh>

Collaboration diagram for G4MaterialScanner:

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)
 
G4String GetRegionName () const
 

Private Member Functions

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

Private Attributes

G4RayShootertheRayShooter
 
G4MatScanMessengertheMessenger
 
G4EventManagertheEventManager
 
G4UserEventActiontheUserEventAction
 
G4UserStackingActiontheUserStackingAction
 
G4UserTrackingActiontheUserTrackingAction
 
G4UserSteppingActiontheUserSteppingAction
 
G4UserEventActiontheMatScannerEventAction
 
G4UserStackingActiontheMatScannerStackingAction
 
G4UserTrackingActiontheMatScannerTrackingAction
 
G4MSSteppingActiontheMatScannerSteppingAction
 
G4ThreeVector eyePosition
 
G4int nTheta
 
G4double thetaMin
 
G4double thetaSpan
 
G4int nPhi
 
G4double phiMin
 
G4double phiSpan
 
G4ThreeVector eyeDirection
 
G4bool regionSensitive
 
G4String regionName
 
G4RegiontheRegion
 

Detailed Description

Definition at line 54 of file G4MaterialScanner.hh.

Constructor & Destructor Documentation

◆ G4MaterialScanner()

G4MaterialScanner::G4MaterialScanner ( )

Definition at line 51 of file G4MaterialScanner.cc.

52 {
54  theMessenger = new G4MatScanMessenger(this);
56 
61 
66 
67  eyePosition = G4ThreeVector(0.,0.,0.);
68  nTheta = 91;
69  thetaMin = 0.*deg;
70  thetaSpan = 90.*deg;
71  nPhi = 37;
72  phiMin = 0.*deg;
73  phiSpan = 360.*deg;
74 
75  regionSensitive = false;
76  regionName = "notDefined";
77  theRegion = 0;
78 }
G4UserEventAction * theMatScannerEventAction
G4UserStackingAction * theUserStackingAction
CLHEP::Hep3Vector G4ThreeVector
G4UserTrackingAction * theMatScannerTrackingAction
G4UserSteppingAction * theUserSteppingAction
static const double deg
Definition: G4SIunits.hh:151
G4UserEventAction * theUserEventAction
G4EventManager * theEventManager
static G4EventManager * GetEventManager()
G4MSSteppingAction * theMatScannerSteppingAction
G4MatScanMessenger * theMessenger
G4UserTrackingAction * theUserTrackingAction
G4UserStackingAction * theMatScannerStackingAction
G4ThreeVector eyePosition
G4RayShooter * theRayShooter
Here is the call graph for this function:

◆ ~G4MaterialScanner()

G4MaterialScanner::~G4MaterialScanner ( )

Definition at line 80 of file G4MaterialScanner.cc.

81 {
82  delete theRayShooter;
84  delete theMessenger;
85 }
G4MSSteppingAction * theMatScannerSteppingAction
G4MatScanMessenger * theMessenger
G4RayShooter * theRayShooter

Member Function Documentation

◆ DoScan()

void G4MaterialScanner::DoScan ( )
private

Definition at line 137 of file G4MaterialScanner.cc.

138 {
139 // Confirm material table is updated
141 
148 
149 // Close geometry and set the application state
151  geomManager->OpenGeometry();
152  geomManager->CloseGeometry(1,0);
153 
154  G4ThreeVector center(0,0,0);
157  navigator->LocateGlobalPointAndSetup(center,0,false);
158 
160  theStateMan->SetNewState(G4State_GeomClosed);
161 
162 // Event loop
163  G4int iEvent = 0;
164  for(G4int iTheta=0;iTheta<nTheta;iTheta++)
165  {
166  G4double theta = thetaMin;
167  if(iTheta>0) theta += G4double(iTheta)*thetaSpan/G4double(nTheta-1);
168  G4double aveLength = 0.;
169  G4double aveX0 = 0.;
170  G4double aveLambda = 0.;
171  G4cout << G4endl;
172  G4cout << " Theta(deg) Phi(deg) Length(mm) x0 lambda0" << G4endl;
173  G4cout << G4endl;
174  for(G4int iPhi=0;iPhi<nPhi;iPhi++)
175  {
176  G4Event* anEvent = new G4Event(iEvent++);
177  G4double phi = phiMin;
178  if(iPhi>0) phi += G4double(iPhi)*phiSpan/G4double(nPhi-1);
179  eyeDirection = G4ThreeVector(std::cos(theta)*std::cos(phi),
180  std::cos(theta)*std::sin(phi),
181  std::sin(theta));
188 
189  G4cout << " "
190  << std::setw(11) << theta/deg << " "
191  << std::setw(11) << phi/deg << " "
192  << std::setw(11) << length/mm << " "
193  << std::setw(11) << x0 << " "
194  << std::setw(11) << lambda << G4endl;
195  aveLength += length/mm;
196  aveX0 += x0;
197  aveLambda += lambda;
198  }
199  if(nPhi>1)
200  {
201  G4cout << G4endl;
202  G4cout << " ave. for theta = " << std::setw(11) << theta/deg << " : "
203  << std::setw(11) << aveLength/nPhi << " "
204  << std::setw(11) << aveX0/nPhi << " "
205  << std::setw(11) << aveLambda/nPhi << G4endl;
206  }
207  }
208 
209  theStateMan->SetNewState(G4State_Idle);
210  return;
211 }
G4ThreeVector eyeDirection
CLHEP::Hep3Vector G4ThreeVector
void Initialize(G4bool rSens, G4Region *reg)
int G4int
Definition: G4Types.hh:78
static G4RunManagerKernel * GetRunManagerKernel()
static G4StateManager * GetStateManager()
G4GLOB_DLL std::ostream G4cout
G4bool SetNewState(G4ApplicationState requestedState)
static const double deg
Definition: G4SIunits.hh:151
G4double GetX0() const
G4Navigator * GetNavigatorForTracking() const
static G4GeometryManager * GetInstance()
static G4TransportationManager * GetTransportationManager()
G4EventManager * theEventManager
void Shoot(G4Event *evt, G4ThreeVector vtx, G4ThreeVector direc)
Definition: G4RayShooter.cc:59
virtual G4VPhysicalVolume * LocateGlobalPointAndSetup(const G4ThreeVector &point, const G4ThreeVector *direction=0, const G4bool pRelativeSearch=true, const G4bool ignoreDirection=true)
Definition: G4Navigator.cc:125
G4MSSteppingAction * theMatScannerSteppingAction
#define G4endl
Definition: G4ios.hh:61
void OpenGeometry(G4VPhysicalVolume *vol=0)
void ProcessOneEvent(G4Event *anEvent)
G4double GetLambda0() const
double G4double
Definition: G4Types.hh:76
G4double GetTotalStepLength() const
static const double mm
Definition: G4SIunits.hh:114
G4ThreeVector eyePosition
G4bool CloseGeometry(G4bool pOptimise=true, G4bool verbose=false, G4VPhysicalVolume *vol=0)
G4RayShooter * theRayShooter
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetEyePosition()

G4ThreeVector G4MaterialScanner::GetEyePosition ( ) const
inline

Definition at line 106 of file G4MaterialScanner.hh.

106 { return eyePosition; }
G4ThreeVector eyePosition
Here is the caller graph for this function:

◆ GetNPhi()

G4int G4MaterialScanner::GetNPhi ( ) const
inline

Definition at line 114 of file G4MaterialScanner.hh.

114 { return nPhi; }
Here is the caller graph for this function:

◆ GetNTheta()

G4int G4MaterialScanner::GetNTheta ( ) const
inline

Definition at line 108 of file G4MaterialScanner.hh.

108 { return nTheta; }
Here is the caller graph for this function:

◆ GetPhiMin()

G4double G4MaterialScanner::GetPhiMin ( ) const
inline

Definition at line 116 of file G4MaterialScanner.hh.

116 { return phiMin; }
Here is the caller graph for this function:

◆ GetPhiSpan()

G4double G4MaterialScanner::GetPhiSpan ( ) const
inline

Definition at line 118 of file G4MaterialScanner.hh.

118 { return phiSpan; }
Here is the caller graph for this function:

◆ GetRegionName()

G4String G4MaterialScanner::GetRegionName ( ) const
inline

Definition at line 122 of file G4MaterialScanner.hh.

122 { return regionName; }
Here is the caller graph for this function:

◆ GetRegionSensitive()

G4bool G4MaterialScanner::GetRegionSensitive ( ) const
inline

Definition at line 120 of file G4MaterialScanner.hh.

120 { return regionSensitive; }
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetThetaMin()

G4double G4MaterialScanner::GetThetaMin ( ) const
inline

Definition at line 110 of file G4MaterialScanner.hh.

110 { return thetaMin; }
Here is the caller graph for this function:

◆ GetThetaSpan()

G4double G4MaterialScanner::GetThetaSpan ( ) const
inline

Definition at line 112 of file G4MaterialScanner.hh.

112 { return thetaSpan; }
Here is the caller graph for this function:

◆ RestoreUserActions()

void G4MaterialScanner::RestoreUserActions ( )
private

Definition at line 125 of file G4MaterialScanner.cc.

126 {
131 
133  if(theSDMan)
134  { theSDMan->Activate("/",true); }
135 }
G4UserStackingAction * theUserStackingAction
G4UserSteppingAction * theUserSteppingAction
void Activate(G4String dName, G4bool activeFlag)
Definition: G4SDManager.cc:121
G4UserEventAction * theUserEventAction
G4EventManager * theEventManager
void SetUserAction(G4UserEventAction *userAction)
G4UserTrackingAction * theUserTrackingAction
static G4SDManager * GetSDMpointerIfExist()
Definition: G4SDManager.cc:49
Here is the call graph for this function:
Here is the caller graph for this function:

◆ Scan()

void G4MaterialScanner::Scan ( )

Definition at line 87 of file G4MaterialScanner.cc.

88 {
90  G4ApplicationState currentState = theStateMan->GetCurrentState();
91  if(currentState!=G4State_Idle)
92  {
93  G4cerr << "Illegal application state - Scan() ignored." << G4endl;
94  return;
95  }
96 
100  DoScan();
102 }
static G4StateManager * GetStateManager()
G4MSSteppingAction * theMatScannerSteppingAction
#define G4endl
Definition: G4ios.hh:61
G4ApplicationState GetCurrentState() const
G4ApplicationState
G4GLOB_DLL std::ostream G4cerr
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetEyePosition()

void G4MaterialScanner::SetEyePosition ( const G4ThreeVector val)
inline

Definition at line 105 of file G4MaterialScanner.hh.

105 { eyePosition = val; }
G4ThreeVector eyePosition
Here is the caller graph for this function:

◆ SetNPhi()

void G4MaterialScanner::SetNPhi ( G4int  val)
inline

Definition at line 113 of file G4MaterialScanner.hh.

113 { nPhi = val; }
Here is the caller graph for this function:

◆ SetNTheta()

void G4MaterialScanner::SetNTheta ( G4int  val)
inline

Definition at line 107 of file G4MaterialScanner.hh.

107 { nTheta = val; }
Here is the caller graph for this function:

◆ SetPhiMin()

void G4MaterialScanner::SetPhiMin ( G4double  val)
inline

Definition at line 115 of file G4MaterialScanner.hh.

115 { phiMin = val; }
Here is the caller graph for this function:

◆ SetPhiSpan()

void G4MaterialScanner::SetPhiSpan ( G4double  val)
inline

Definition at line 117 of file G4MaterialScanner.hh.

117 { phiSpan = val; }
Here is the caller graph for this function:

◆ SetRegionName()

G4bool G4MaterialScanner::SetRegionName ( const G4String val)

Definition at line 213 of file G4MaterialScanner.cc.

214 {
215  G4Region* aRegion = G4RegionStore::GetInstance()->GetRegion(val);
216  if(aRegion)
217  {
218  theRegion = aRegion;
219  regionName = val;
220  return true;
221  }
222  else
223  {
224  G4cerr << "Region <" << val << "> not found. Command ignored." << G4endl;
225  G4cerr << "Defined regions are : " << G4endl;
226  for(size_t i=0;i<G4RegionStore::GetInstance()->size();i++)
227  { G4cerr << " " << (*(G4RegionStore::GetInstance()))[i]->GetName(); }
228  G4cerr << G4endl;
229  return false;
230  }
231 }
static G4RegionStore * GetInstance()
#define G4endl
Definition: G4ios.hh:61
G4Region * GetRegion(const G4String &name, G4bool verbose=true) const
G4GLOB_DLL std::ostream G4cerr
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetRegionSensitive()

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

Definition at line 119 of file G4MaterialScanner.hh.

119 { regionSensitive = val; }
Here is the caller graph for this function:

◆ SetThetaMin()

void G4MaterialScanner::SetThetaMin ( G4double  val)
inline

Definition at line 109 of file G4MaterialScanner.hh.

109 { thetaMin = val; }
Here is the caller graph for this function:

◆ SetThetaSpan()

void G4MaterialScanner::SetThetaSpan ( G4double  val)
inline

Definition at line 111 of file G4MaterialScanner.hh.

111 { thetaSpan = val; }
Here is the caller graph for this function:

◆ StoreUserActions()

void G4MaterialScanner::StoreUserActions ( )
private

Definition at line 104 of file G4MaterialScanner.cc.

105 {
110 
115 
117  if(theSDMan)
118  { theSDMan->Activate("/",false); }
119 
121  theGeomMan->OpenGeometry();
122  theGeomMan->CloseGeometry(true);
123 }
G4UserEventAction * theMatScannerEventAction
G4UserStackingAction * theUserStackingAction
G4UserTrackingAction * theMatScannerTrackingAction
G4UserEventAction * GetUserEventAction()
G4UserSteppingAction * theUserSteppingAction
void Activate(G4String dName, G4bool activeFlag)
Definition: G4SDManager.cc:121
G4UserSteppingAction * GetUserSteppingAction()
static G4GeometryManager * GetInstance()
G4UserEventAction * theUserEventAction
G4EventManager * theEventManager
void SetUserAction(G4UserEventAction *userAction)
G4UserTrackingAction * GetUserTrackingAction()
G4MSSteppingAction * theMatScannerSteppingAction
void OpenGeometry(G4VPhysicalVolume *vol=0)
G4UserTrackingAction * theUserTrackingAction
G4UserStackingAction * theMatScannerStackingAction
G4bool CloseGeometry(G4bool pOptimise=true, G4bool verbose=false, G4VPhysicalVolume *vol=0)
G4UserStackingAction * GetUserStackingAction()
static G4SDManager * GetSDMpointerIfExist()
Definition: G4SDManager.cc:49
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ eyeDirection

G4ThreeVector G4MaterialScanner::eyeDirection
private

Definition at line 98 of file G4MaterialScanner.hh.

◆ eyePosition

G4ThreeVector G4MaterialScanner::eyePosition
private

Definition at line 90 of file G4MaterialScanner.hh.

◆ nPhi

G4int G4MaterialScanner::nPhi
private

Definition at line 94 of file G4MaterialScanner.hh.

◆ nTheta

G4int G4MaterialScanner::nTheta
private

Definition at line 91 of file G4MaterialScanner.hh.

◆ phiMin

G4double G4MaterialScanner::phiMin
private

Definition at line 95 of file G4MaterialScanner.hh.

◆ phiSpan

G4double G4MaterialScanner::phiSpan
private

Definition at line 96 of file G4MaterialScanner.hh.

◆ regionName

G4String G4MaterialScanner::regionName
private

Definition at line 101 of file G4MaterialScanner.hh.

◆ regionSensitive

G4bool G4MaterialScanner::regionSensitive
private

Definition at line 100 of file G4MaterialScanner.hh.

◆ theEventManager

G4EventManager* G4MaterialScanner::theEventManager
private

Definition at line 78 of file G4MaterialScanner.hh.

◆ theMatScannerEventAction

G4UserEventAction* G4MaterialScanner::theMatScannerEventAction
private

Definition at line 85 of file G4MaterialScanner.hh.

◆ theMatScannerStackingAction

G4UserStackingAction* G4MaterialScanner::theMatScannerStackingAction
private

Definition at line 86 of file G4MaterialScanner.hh.

◆ theMatScannerSteppingAction

G4MSSteppingAction* G4MaterialScanner::theMatScannerSteppingAction
private

Definition at line 88 of file G4MaterialScanner.hh.

◆ theMatScannerTrackingAction

G4UserTrackingAction* G4MaterialScanner::theMatScannerTrackingAction
private

Definition at line 87 of file G4MaterialScanner.hh.

◆ theMessenger

G4MatScanMessenger* G4MaterialScanner::theMessenger
private

Definition at line 76 of file G4MaterialScanner.hh.

◆ theRayShooter

G4RayShooter* G4MaterialScanner::theRayShooter
private

Definition at line 75 of file G4MaterialScanner.hh.

◆ theRegion

G4Region* G4MaterialScanner::theRegion
private

Definition at line 102 of file G4MaterialScanner.hh.

◆ thetaMin

G4double G4MaterialScanner::thetaMin
private

Definition at line 92 of file G4MaterialScanner.hh.

◆ thetaSpan

G4double G4MaterialScanner::thetaSpan
private

Definition at line 93 of file G4MaterialScanner.hh.

◆ theUserEventAction

G4UserEventAction* G4MaterialScanner::theUserEventAction
private

Definition at line 80 of file G4MaterialScanner.hh.

◆ theUserStackingAction

G4UserStackingAction* G4MaterialScanner::theUserStackingAction
private

Definition at line 81 of file G4MaterialScanner.hh.

◆ theUserSteppingAction

G4UserSteppingAction* G4MaterialScanner::theUserSteppingAction
private

Definition at line 83 of file G4MaterialScanner.hh.

◆ theUserTrackingAction

G4UserTrackingAction* G4MaterialScanner::theUserTrackingAction
private

Definition at line 82 of file G4MaterialScanner.hh.


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