Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
WLSSteppingAction.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
28 //
29 //
30 //
31 
32 #include "G4Run.hh"
33 #include "G4Step.hh"
34 #include "G4Track.hh"
35 #include "G4StepPoint.hh"
36 #include "G4TrackStatus.hh"
37 #include "G4VPhysicalVolume.hh"
38 #include "G4ParticleDefinition.hh"
39 
40 #include "WLSSteppingAction.hh"
43 #include "WLSPhotonDetSD.hh"
44 
45 #include "G4ParticleTypes.hh"
46 
48 
49 #include "G4ProcessManager.hh"
50 #include "G4OpBoundaryProcess.hh"
51 
52 #include "G4RunManager.hh"
53 #include "G4SDManager.hh"
54 #include "G4UImanager.hh"
55 
56 #include "G4ThreeVector.hh"
57 #include "G4ios.hh"
58 #include "G4SystemOfUnits.hh"
59 #include <sstream>
60 
61 // Purpose: Save relevant information into User Track Information
62 
63 static const G4ThreeVector ZHat = G4ThreeVector(0.0,0.0,1.0);
64 
65 G4int WLSSteppingAction::maxRndmSave = 10000;
66 
68  : detector(DC)
69 {
70  steppingMessenger = new WLSSteppingActionMessenger(this);
71 
72  counterEnd = 0;
73  counterMid = 0;
74  bounceLimit = 100000;
75 
76  ResetCounters();
77 }
78 
80 {
81  delete steppingMessenger;
82 }
83 
84 void WLSSteppingAction::SetBounceLimit(G4int i) {bounceLimit = i;}
86 G4int WLSSteppingAction::GetNumberOfClad1Bounces() {return counterClad1Bounce;}
87 G4int WLSSteppingAction::GetNumberOfClad2Bounces() {return counterClad2Bounce;}
88 G4int WLSSteppingAction::GetNumberOfWLSBounces() {return counterWLSBounce;}
90  G4int temp = counterEnd; counterEnd = 0; return temp;
91 }
92 
93 // save the random status into a sub-directory
94 // Pre: subDir must be empty or ended with "/"
95 inline void WLSSteppingAction::saveRandomStatus(G4String subDir) {
96 
97  // don't save if the maximum amount has been reached
98  if (WLSSteppingAction::maxRndmSave == 0) return;
99 
100  G4RunManager* theRunManager = G4RunManager::GetRunManager();
101  G4String randomNumberStatusDir = theRunManager->GetRandomNumberStoreDir();
102 
103  G4String fileIn = randomNumberStatusDir + "currentEvent.rndm";
104 
105  std::ostringstream os;
106 
107  os << "run" << theRunManager->GetCurrentRun()->GetRunID() << "evt"
108  << theRunManager->GetCurrentEvent()->GetEventID() << ".rndm" << '\0';
109 
110  G4String fileOut = randomNumberStatusDir + subDir + os.str();
111 
112  G4String copCmd = "/control/shell cp "+fileIn+" "+fileOut;
114 
115  WLSSteppingAction::maxRndmSave--;
116 }
117 
118 void WLSSteppingAction::UpdateHistogramSuccess(G4StepPoint* ,G4Track* ) {}
119 
120 void WLSSteppingAction::UpdateHistogramReflect(G4StepPoint* ,G4Track* ) {}
121 
122 void WLSSteppingAction::UpdateHistogramEscape(G4StepPoint* , G4Track* ) {}
123 
124 void WLSSteppingAction::UpdateHistogramAbsorb(G4StepPoint* , G4Track* ) {}
125 
127 {
128  G4Track* theTrack = theStep->GetTrack();
129  WLSUserTrackInformation* trackInformation
131 
132  G4StepPoint* thePrePoint = theStep->GetPreStepPoint();
133  G4StepPoint* thePostPoint = theStep->GetPostStepPoint();
134 
135  G4VPhysicalVolume* thePrePV = thePrePoint->GetPhysicalVolume();
136  G4VPhysicalVolume* thePostPV = thePostPoint->GetPhysicalVolume();
137 
138  G4String thePrePVname = " ";
139  G4String thePostPVname = " ";
140 
141  if (thePostPV) {
142  thePrePVname = thePrePV->GetName();
143  thePostPVname = thePostPV->GetName();
144  }
145 
146  //Recording data for start
147  if (theTrack->GetParentID()==0) {
148  //This is a primary track
149  if ( theTrack->GetCurrentStepNumber() == 1 ) {
150 // G4double x = theTrack->GetVertexPosition().x();
151 // G4double y = theTrack->GetVertexPosition().y();
152  G4double z = theTrack->GetVertexPosition().z();
153 // G4double pz = theTrack->GetVertexMomentumDirection().z();
154  initTheta = theTrack->GetVertexMomentumDirection().angle(ZHat);
155  initZ = z;
156  }
157  }
158 
159  // Retrieve the status of the photon
161 
162  G4ProcessManager* OpManager =
164 
165  if (OpManager) {
166  G4int MAXofPostStepLoops =
167  OpManager->GetPostStepProcessVector()->entries();
168  G4ProcessVector* fPostStepDoItVector =
170 
171  for ( G4int i=0; i<MAXofPostStepLoops; i++) {
172  G4VProcess* fCurrentProcess = (*fPostStepDoItVector)[i];
173  opProcess = dynamic_cast<G4OpBoundaryProcess*>(fCurrentProcess);
174  if (opProcess) { theStatus = opProcess->GetStatus(); break;}
175  }
176  }
177 
178  // Find the skewness of the ray at first change of boundary
179  if ( initGamma == -1 &&
180  (theStatus == TotalInternalReflection
181  || theStatus == FresnelReflection
182  || theStatus == FresnelRefraction)
183  && trackInformation->isStatus(InsideOfFiber) ) {
184 
185  G4double px = theTrack->GetVertexMomentumDirection().x();
186  G4double py = theTrack->GetVertexMomentumDirection().y();
187  G4double x = theTrack->GetPosition().x();
188  G4double y = theTrack->GetPosition().y();
189 
190  initGamma = x * px + y * py;
191 
192  initGamma = initGamma / std::sqrt(px*px + py*py) / std::sqrt(x*x + y*y);
193 
194  initGamma = std::acos(initGamma*rad);
195 
196  if ( initGamma / deg > 90.0) { initGamma = 180 * deg - initGamma;}
197  }
198  // Record Photons that missed the photon detector but escaped from readout
199  if ( !thePostPV && trackInformation->isStatus(EscapedFromReadOut) ) {
200 // UpdateHistogramSuccess(thePostPoint,theTrack);
201  ResetCounters();
202 
203  return;
204  }
205 
206  // Assumed photons are originated at the fiber OR
207  // the fiber is the first material the photon hits
208  switch (theStatus) {
209 
210  // Exiting the fiber
211  case FresnelRefraction:
212  case SameMaterial:
213 
214  G4bool isFiber;
215  isFiber = thePostPVname == "WLSFiber"
216  || thePostPVname == "Clad1"
217  || thePostPVname == "Clad2";
218 
219  if ( isFiber ) {
220 
221  if (trackInformation->isStatus(OutsideOfFiber))
222  trackInformation->AddStatusFlag(InsideOfFiber);
223 
224  // Set the Exit flag when the photon refracted out of the fiber
225  } else if (trackInformation->isStatus(InsideOfFiber)) {
226 
227  // EscapedFromReadOut if the z position is the same as fiber's end
228  if (theTrack->GetPosition().z() == detector->GetWLSFiberEnd())
229  {
230  trackInformation->AddStatusFlag(EscapedFromReadOut);
231  counterEnd++;
232  }
233  else // Escaped from side
234  {
235  trackInformation->AddStatusFlag(EscapedFromSide);
236  trackInformation->SetExitPosition(theTrack->GetPosition());
237 
238 // UpdateHistogramEscape(thePostPoint,theTrack);
239 
240  counterMid++;
241  ResetCounters();
242  }
243 
244  trackInformation->AddStatusFlag(OutsideOfFiber);
245  trackInformation->SetExitPosition(theTrack->GetPosition());
246 
247  }
248 
249  return;
250 
251  // Internal Reflections
253 
254  // Kill the track if it's number of bounces exceeded the limit
255  if (bounceLimit > 0 && counterBounce >= bounceLimit)
256  {
257  theTrack->SetTrackStatus(fStopAndKill);
258  trackInformation->AddStatusFlag(murderee);
259  ResetCounters();
260  G4cout << "\n Bounce Limit Exceeded" << G4endl;
261  return;
262  }
263 
264  case FresnelReflection:
265 
266  counterBounce++;
267 
268  if ( thePrePVname == "WLSFiber") counterWLSBounce++;
269 
270  else if ( thePrePVname == "Clad1") counterClad1Bounce++;
271 
272  else if ( thePrePVname == "Clad2") counterClad2Bounce++;
273 
274  // Determine if the photon has reflected off the read-out end
275  if (theTrack->GetPosition().z() == detector->GetWLSFiberEnd())
276  {
277  if (!trackInformation->isStatus(ReflectedAtReadOut) &&
278  trackInformation->isStatus(InsideOfFiber))
279  {
280  trackInformation->AddStatusFlag(ReflectedAtReadOut);
281 
282  if (detector->IsPerfectFiber() &&
283  theStatus == TotalInternalReflection)
284  {
285  theTrack->SetTrackStatus(fStopAndKill);
286  trackInformation->AddStatusFlag(murderee);
287 // UpdateHistogramReflect(thePostPoint,theTrack);
288  ResetCounters();
289  return;
290  }
291  }
292  }
293  return;
294 
295  // Reflection of the mirror
297  case LobeReflection:
298  case SpikeReflection:
299 
300  // Check if it hits the mirror
301  if ( thePostPVname == "Mirror" )
302  trackInformation->AddStatusFlag(ReflectedAtMirror);
303 
304  return;
305 
306  // Detected by a detector
307  case Detection:
308 
309  // Check if the photon hits the detector and process the hit if it does
310  if ( thePostPVname == "PhotonDet" ) {
311 
313  G4String SDname="WLS/PhotonDet";
314  WLSPhotonDetSD* mppcSD =
315  (WLSPhotonDetSD*)SDman->FindSensitiveDetector(SDname);
316 
317  if (mppcSD) mppcSD->ProcessHits_constStep(theStep,NULL);
318 
319  // Record Photons that escaped at the end
320 // if (trackInformation->isStatus(EscapedFromReadOut))
321 // UpdateHistogramSuccess(thePostPoint,theTrack);
322 
323  // Stop Tracking when it hits the detector's surface
324  ResetCounters();
325  theTrack->SetTrackStatus(fStopAndKill);
326 
327  return;
328  }
329 
330  break;
331 
332  default: break;
333 
334  }
335 
336  // Check for absorbed photons
337  if (theTrack->GetTrackStatus() != fAlive &&
338  trackInformation->isStatus(InsideOfFiber))
339  {
340 // UpdateHistogramAbsorb(thePostPoint,theTrack);
341  ResetCounters();
342  return;
343  }
344 
345 }