Geant4  10.02.p03
WLSSteppingAction Class Reference

#include <WLSSteppingAction.hh>

Inheritance diagram for WLSSteppingAction:
Collaboration diagram for WLSSteppingAction:

Public Member Functions

 WLSSteppingAction (WLSDetectorConstruction *)
 
virtual ~WLSSteppingAction ()
 
virtual void UserSteppingAction (const G4Step *)
 
void SetBounceLimit (G4int)
 
G4int GetNumberOfBounces ()
 
G4int GetNumberOfClad1Bounces ()
 
G4int GetNumberOfClad2Bounces ()
 
G4int GetNumberOfWLSBounces ()
 
G4int ResetSuccessCounter ()
 
- Public Member Functions inherited from G4UserSteppingAction
 G4UserSteppingAction ()
 
virtual ~G4UserSteppingAction ()
 
void SetSteppingManagerPointer (G4SteppingManager *pValue)
 

Private Member Functions

void ResetCounters ()
 
void saveRandomStatus (G4String subDir)
 

Private Attributes

G4int fBounceLimit
 
G4int fCounterEnd
 
G4int fCounterMid
 
G4int fCounterBounce
 
G4int fCounterWLSBounce
 
G4int fCounterClad1Bounce
 
G4int fCounterClad2Bounce
 
G4double fInitGamma
 
G4double fInitTheta
 
G4OpBoundaryProcessfOpProcess
 
WLSDetectorConstructionfDetector
 
WLSSteppingActionMessengerfSteppingMessenger
 

Static Private Attributes

static G4int fMaxRndmSave = 10000
 

Additional Inherited Members

- Protected Attributes inherited from G4UserSteppingAction
G4SteppingManagerfpSteppingManager
 

Detailed Description

Definition at line 48 of file WLSSteppingAction.hh.

Constructor & Destructor Documentation

◆ WLSSteppingAction()

WLSSteppingAction::WLSSteppingAction ( WLSDetectorConstruction detector)

Definition at line 69 of file WLSSteppingAction.cc.

70  : fDetector(detector)
71 {
73 
74  fCounterEnd = 0;
75  fCounterMid = 0;
76  fBounceLimit = 100000;
77 
78  fOpProcess = NULL;
79 
80  ResetCounters();
81 }
G4OpBoundaryProcess * fOpProcess
WLSSteppingActionMessenger * fSteppingMessenger
WLSDetectorConstruction * fDetector
Here is the call graph for this function:

◆ ~WLSSteppingAction()

WLSSteppingAction::~WLSSteppingAction ( )
virtual

Definition at line 85 of file WLSSteppingAction.cc.

86 {
87  delete fSteppingMessenger;
88 }
WLSSteppingActionMessenger * fSteppingMessenger

Member Function Documentation

◆ GetNumberOfBounces()

G4int WLSSteppingAction::GetNumberOfBounces ( )

Definition at line 96 of file WLSSteppingAction.cc.

96 {return fCounterBounce;}

◆ GetNumberOfClad1Bounces()

G4int WLSSteppingAction::GetNumberOfClad1Bounces ( )

Definition at line 100 of file WLSSteppingAction.cc.

100 {return fCounterClad1Bounce;}

◆ GetNumberOfClad2Bounces()

G4int WLSSteppingAction::GetNumberOfClad2Bounces ( )

Definition at line 104 of file WLSSteppingAction.cc.

104 {return fCounterClad2Bounce;}

◆ GetNumberOfWLSBounces()

G4int WLSSteppingAction::GetNumberOfWLSBounces ( )

Definition at line 108 of file WLSSteppingAction.cc.

108 {return fCounterWLSBounce;}

◆ ResetCounters()

void WLSSteppingAction::ResetCounters ( )
inlineprivate

Definition at line 98 of file WLSSteppingAction.hh.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ ResetSuccessCounter()

G4int WLSSteppingAction::ResetSuccessCounter ( )

Definition at line 112 of file WLSSteppingAction.cc.

112  {
113  G4int temp = fCounterEnd; fCounterEnd = 0; return temp;
114 }
int G4int
Definition: G4Types.hh:78

◆ saveRandomStatus()

void WLSSteppingAction::saveRandomStatus ( G4String  subDir)
inlineprivate

Definition at line 118 of file WLSSteppingAction.cc.

121 {
122 
123  // don't save if the maximum amount has been reached
124  if (WLSSteppingAction::fMaxRndmSave == 0) return;
125 
126  G4RunManager* theRunManager = G4RunManager::GetRunManager();
127  G4String randomNumberStatusDir = theRunManager->GetRandomNumberStoreDir();
128 
129  G4String fileIn = randomNumberStatusDir + "currentEvent.rndm";
130 
131  std::ostringstream os;
132 
133  os << "run" << theRunManager->GetCurrentRun()->GetRunID() << "evt"
134  << theRunManager->GetCurrentEvent()->GetEventID() << ".rndm" << '\0';
135 
136  G4String fileOut = randomNumberStatusDir + subDir + os.str();
137 
138  G4String copCmd = "/control/shell cp "+fileIn+" "+fileOut;
140 
142 }
const G4String & GetRandomNumberStoreDir() const
static G4UImanager * GetUIpointer()
Definition: G4UImanager.cc:58
static G4int fMaxRndmSave
const G4Run * GetCurrentRun() const
G4int GetRunID() const
Definition: G4Run.hh:76
const G4Event * GetCurrentEvent() const
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:79
G4int GetEventID() const
Definition: G4Event.hh:151
G4int ApplyCommand(const char *aCommand)
Definition: G4UImanager.cc:446
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetBounceLimit()

void WLSSteppingAction::SetBounceLimit ( G4int  i)

Definition at line 92 of file WLSSteppingAction.cc.

◆ UserSteppingAction()

void WLSSteppingAction::UserSteppingAction ( const G4Step *  theStep)
virtual

Reimplemented from G4UserSteppingAction.

Definition at line 146 of file WLSSteppingAction.cc.

147 {
148  G4Track* theTrack = theStep->GetTrack();
149  WLSUserTrackInformation* trackInformation
150  = (WLSUserTrackInformation*)theTrack->GetUserInformation();
151 
152  G4StepPoint* thePrePoint = theStep->GetPreStepPoint();
153  G4StepPoint* thePostPoint = theStep->GetPostStepPoint();
154 
155  G4VPhysicalVolume* thePrePV = thePrePoint->GetPhysicalVolume();
156  G4VPhysicalVolume* thePostPV = thePostPoint->GetPhysicalVolume();
157 
158  G4String thePrePVname = " ";
159  G4String thePostPVname = " ";
160 
161  if (thePostPV) {
162  thePrePVname = thePrePV->GetName();
163  thePostPVname = thePostPV->GetName();
164  }
165 
166  //Recording data for start
167  if (theTrack->GetParentID()==0) {
168  //This is a primary track
169  if ( theTrack->GetCurrentStepNumber() == 1 ) {
170 // G4double x = theTrack->GetVertexPosition().x();
171 // G4double y = theTrack->GetVertexPosition().y();
172 // G4double z = theTrack->GetVertexPosition().z();
173 // G4double pz = theTrack->GetVertexMomentumDirection().z();
174 // G4double fInitTheta =
175 // theTrack->GetVertexMomentumDirection().angle(ZHat);
176  }
177  }
178 
179  // Retrieve the status of the photon
181 
182  G4ProcessManager* OpManager =
184 
185  if (OpManager) {
186  G4int MAXofPostStepLoops =
187  OpManager->GetPostStepProcessVector()->entries();
188  G4ProcessVector* fPostStepDoItVector =
190 
191  for ( G4int i=0; i<MAXofPostStepLoops; i++) {
192  G4VProcess* fCurrentProcess = (*fPostStepDoItVector)[i];
193  fOpProcess = dynamic_cast<G4OpBoundaryProcess*>(fCurrentProcess);
194  if (fOpProcess) { theStatus = fOpProcess->GetStatus(); break;}
195  }
196  }
197 
198  // Find the skewness of the ray at first change of boundary
199  if ( fInitGamma == -1 &&
200  (theStatus == TotalInternalReflection
201  || theStatus == FresnelReflection
202  || theStatus == FresnelRefraction)
203  && trackInformation->isStatus(InsideOfFiber) ) {
204 
205  G4double px = theTrack->GetVertexMomentumDirection().x();
206  G4double py = theTrack->GetVertexMomentumDirection().y();
207  G4double x = theTrack->GetPosition().x();
208  G4double y = theTrack->GetPosition().y();
209 
210  fInitGamma = x * px + y * py;
211 
212  fInitGamma =
213  fInitGamma / std::sqrt(px*px + py*py) / std::sqrt(x*x + y*y);
214 
215  fInitGamma = std::acos(fInitGamma*rad);
216 
217  if ( fInitGamma / deg > 90.0) { fInitGamma = 180 * deg - fInitGamma;}
218  }
219  // Record Photons that missed the photon detector but escaped from readout
220  if ( !thePostPV && trackInformation->isStatus(EscapedFromReadOut) ) {
221 // UpdateHistogramSuccess(thePostPoint,theTrack);
222  ResetCounters();
223 
224  return;
225  }
226 
227  // Assumed photons are originated at the fiber OR
228  // the fiber is the first material the photon hits
229  switch (theStatus) {
230 
231  // Exiting the fiber
232  case FresnelRefraction:
233  case SameMaterial:
234 
235  G4bool isFiber;
236  isFiber = thePostPVname == "WLSFiber"
237  || thePostPVname == "Clad1"
238  || thePostPVname == "Clad2";
239 
240  if ( isFiber ) {
241 
242  if (trackInformation->isStatus(OutsideOfFiber))
243  trackInformation->AddStatusFlag(InsideOfFiber);
244 
245  // Set the Exit flag when the photon refracted out of the fiber
246  } else if (trackInformation->isStatus(InsideOfFiber)) {
247 
248  // EscapedFromReadOut if the z position is the same as fiber's end
249  if (theTrack->GetPosition().z() == fDetector->GetWLSFiberEnd())
250  {
251  trackInformation->AddStatusFlag(EscapedFromReadOut);
252  fCounterEnd++;
253  }
254  else // Escaped from side
255  {
256  trackInformation->AddStatusFlag(EscapedFromSide);
257  trackInformation->SetExitPosition(theTrack->GetPosition());
258 
259 // UpdateHistogramEscape(thePostPoint,theTrack);
260 
261  fCounterMid++;
262  ResetCounters();
263  }
264 
265  trackInformation->AddStatusFlag(OutsideOfFiber);
266  trackInformation->SetExitPosition(theTrack->GetPosition());
267 
268  }
269 
270  return;
271 
272  // Internal Reflections
274 
275  // Kill the track if it's number of bounces exceeded the limit
277  {
278  theTrack->SetTrackStatus(fStopAndKill);
279  trackInformation->AddStatusFlag(murderee);
280  ResetCounters();
281  G4cout << "\n Bounce Limit Exceeded" << G4endl;
282  return;
283  }
284 
285  case FresnelReflection:
286 
287  fCounterBounce++;
288 
289  if ( thePrePVname == "WLSFiber") fCounterWLSBounce++;
290 
291  else if ( thePrePVname == "Clad1") fCounterClad1Bounce++;
292 
293  else if ( thePrePVname == "Clad2") fCounterClad2Bounce++;
294 
295  // Determine if the photon has reflected off the read-out end
296  if (theTrack->GetPosition().z() == fDetector->GetWLSFiberEnd())
297  {
298  if (!trackInformation->isStatus(ReflectedAtReadOut) &&
299  trackInformation->isStatus(InsideOfFiber))
300  {
301  trackInformation->AddStatusFlag(ReflectedAtReadOut);
302 
303  if (fDetector->IsPerfectFiber() &&
304  theStatus == TotalInternalReflection)
305  {
306  theTrack->SetTrackStatus(fStopAndKill);
307  trackInformation->AddStatusFlag(murderee);
308 // UpdateHistogramReflect(thePostPoint,theTrack);
309  ResetCounters();
310  return;
311  }
312  }
313  }
314  return;
315 
316  // Reflection of the mirror
318  case LobeReflection:
319  case SpikeReflection:
320 
321  // Check if it hits the mirror
322  if ( thePostPVname == "Mirror" )
323  trackInformation->AddStatusFlag(ReflectedAtMirror);
324 
325  return;
326 
327  // Detected by a detector
328  case Detection:
329 
330  // Check if the photon hits the detector and process the hit if it does
331  if ( thePostPVname == "PhotonDet" ) {
332 
334  G4String SDname="WLS/PhotonDet";
335  WLSPhotonDetSD* mppcSD =
336  (WLSPhotonDetSD*)SDman->FindSensitiveDetector(SDname);
337 
338  if (mppcSD) mppcSD->ProcessHits_constStep(theStep,NULL);
339 
340  // Record Photons that escaped at the end
341 // if (trackInformation->isStatus(EscapedFromReadOut))
342 // UpdateHistogramSuccess(thePostPoint,theTrack);
343 
344  // Stop Tracking when it hits the detector's surface
345  ResetCounters();
346  theTrack->SetTrackStatus(fStopAndKill);
347 
348  return;
349  }
350 
351  break;
352 
353  default: break;
354 
355  }
356 
357  // Check for absorbed photons
358  if (theTrack->GetTrackStatus() != fAlive &&
359  trackInformation->isStatus(InsideOfFiber))
360  {
361 // UpdateHistogramAbsorb(thePostPoint,theTrack);
362  ResetCounters();
363  return;
364  }
365 
366 }
G4bool AddStatusFlag(TrackStatus s)
G4bool ProcessHits_constStep(const G4Step *, G4TouchableHistory *)
void SetExitPosition(const G4ThreeVector &pos)
G4OpBoundaryProcessStatus
G4ProcessManager * GetProcessManager() const
G4ProcessVector * GetPostStepProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
int G4int
Definition: G4Types.hh:78
G4int entries() const
Double_t y
G4GLOB_DLL std::ostream G4cout
static const double deg
Definition: G4SIunits.hh:151
bool G4bool
Definition: G4Types.hh:79
G4OpBoundaryProcess * fOpProcess
G4VSensitiveDetector * FindSensitiveDetector(G4String dName, G4bool warning=true)
Definition: G4SDManager.cc:128
static const double rad
Definition: G4SIunits.hh:148
G4bool isStatus(TrackStatus s)
static G4OpticalPhoton * OpticalPhoton()
const G4String & GetName() const
static G4SDManager * GetSDMpointer()
Definition: G4SDManager.cc:40
G4OpBoundaryProcessStatus GetStatus() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
WLSDetectorConstruction * fDetector
Here is the call graph for this function:

Member Data Documentation

◆ fBounceLimit

G4int WLSSteppingAction::fBounceLimit
private

Definition at line 70 of file WLSSteppingAction.hh.

◆ fCounterBounce

G4int WLSSteppingAction::fCounterBounce
private

Definition at line 76 of file WLSSteppingAction.hh.

◆ fCounterClad1Bounce

G4int WLSSteppingAction::fCounterClad1Bounce
private

Definition at line 80 of file WLSSteppingAction.hh.

◆ fCounterClad2Bounce

G4int WLSSteppingAction::fCounterClad2Bounce
private

Definition at line 82 of file WLSSteppingAction.hh.

◆ fCounterEnd

G4int WLSSteppingAction::fCounterEnd
private

Definition at line 72 of file WLSSteppingAction.hh.

◆ fCounterMid

G4int WLSSteppingAction::fCounterMid
private

Definition at line 74 of file WLSSteppingAction.hh.

◆ fCounterWLSBounce

G4int WLSSteppingAction::fCounterWLSBounce
private

Definition at line 78 of file WLSSteppingAction.hh.

◆ fDetector

WLSDetectorConstruction* WLSSteppingAction::fDetector
private

Definition at line 94 of file WLSSteppingAction.hh.

◆ fInitGamma

G4double WLSSteppingAction::fInitGamma
private

Definition at line 85 of file WLSSteppingAction.hh.

◆ fInitTheta

G4double WLSSteppingAction::fInitTheta
private

Definition at line 87 of file WLSSteppingAction.hh.

◆ fMaxRndmSave

G4int WLSSteppingAction::fMaxRndmSave = 10000
staticprivate

Definition at line 92 of file WLSSteppingAction.hh.

◆ fOpProcess

G4OpBoundaryProcess* WLSSteppingAction::fOpProcess
private

Definition at line 89 of file WLSSteppingAction.hh.

◆ fSteppingMessenger

WLSSteppingActionMessenger* WLSSteppingAction::fSteppingMessenger
private

Definition at line 96 of file WLSSteppingAction.hh.


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