Loading...
Searching...
No Matches
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//
26//
27/// \file optical/wls/src/WLSSteppingAction.cc
28/// \brief Implementation of the WLSSteppingAction class
29//
30//
31#include "WLSSteppingAction.hh"
32
34#include "WLSEventAction.hh"
35#include "WLSPhotonDetSD.hh"
38
39#include "G4ios.hh"
40#include "G4OpBoundaryProcess.hh"
41#include "G4OpticalPhoton.hh"
42#include "G4ProcessManager.hh"
43#include "G4Run.hh"
44#include "G4SDManager.hh"
45#include "G4Step.hh"
46#include "G4StepPoint.hh"
47#include "G4SystemOfUnits.hh"
48#include "G4ThreeVector.hh"
49#include "G4Track.hh"
50#include "G4TrackStatus.hh"
51#include "G4UImanager.hh"
52#include "G4VPhysicalVolume.hh"
53
54// Purpose: Save relevant information into User Track Information
55
56//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
57
59 WLSEventAction* event)
60 : fDetector(detector)
61 , fEventAction(event)
62{
65}
66
67//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
68
70
71//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
72
74
75//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
76
78
79//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
80
85
86//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
87
92
93//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
94
96
97//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
98
100{
101 G4int temp = fCounterEnd;
102 fCounterEnd = 0;
103 return temp;
104}
105
106//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
107
109{
110 G4Track* theTrack = theStep->GetTrack();
111 auto trackInformation =
112 (WLSUserTrackInformation*) theTrack->GetUserInformation();
113
114 G4StepPoint* thePrePoint = theStep->GetPreStepPoint();
115 G4StepPoint* thePostPoint = theStep->GetPostStepPoint();
116
117 G4VPhysicalVolume* thePrePV = thePrePoint->GetPhysicalVolume();
118 G4VPhysicalVolume* thePostPV = thePostPoint->GetPhysicalVolume();
119
120 G4String thePrePVname = " ";
121 G4String thePostPVname = " ";
122
123 if(thePostPV)
124 {
125 thePrePVname = thePrePV->GetName();
126 thePostPVname = thePostPV->GetName();
127 }
128
129 // Recording data for start
130 // static const G4ThreeVector ZHat = G4ThreeVector(0.0,0.0,1.0);
131 if(theTrack->GetParentID() == 0)
132 {
133 // This is a primary track
134 if(theTrack->GetCurrentStepNumber() == 1)
135 {
136 // G4double x = theTrack->GetVertexPosition().x();
137 // G4double y = theTrack->GetVertexPosition().y();
138 // G4double z = theTrack->GetVertexPosition().z();
139 // G4double pz = theTrack->GetVertexMomentumDirection().z();
140 // G4double fInitTheta =
141 // theTrack->GetVertexMomentumDirection().angle(ZHat);
142 }
143 }
144
145 // Retrieve the status of the photon
146 G4OpBoundaryProcessStatus theStatus = Undefined;
147
148 static G4ThreadLocal G4ProcessManager* OpManager =
149 G4OpticalPhoton::OpticalPhoton()->GetProcessManager();
150
151 if(OpManager)
152 {
153 G4int nproc = OpManager->GetPostStepProcessVector()->entries();
154 G4ProcessVector* fPostStepDoItVector =
155 OpManager->GetPostStepProcessVector(typeDoIt);
156
157 for(G4int i = 0; i < nproc; ++i)
158 {
159 G4VProcess* fCurrentProcess = (*fPostStepDoItVector)[i];
160 fOpProcess = dynamic_cast<G4OpBoundaryProcess*>(fCurrentProcess);
161 if(fOpProcess)
162 {
163 theStatus = fOpProcess->GetStatus();
164 break;
165 }
166 }
167 }
168
169 // Find the skewness of the ray at first change of boundary
170 if(fInitGamma == -1 &&
171 (theStatus == TotalInternalReflection || theStatus == FresnelReflection ||
172 theStatus == FresnelRefraction) &&
173 trackInformation->IsStatus(InsideOfFiber))
174 {
175 G4double px = theTrack->GetVertexMomentumDirection().x();
176 G4double py = theTrack->GetVertexMomentumDirection().y();
177 G4double x = theTrack->GetPosition().x();
178 G4double y = theTrack->GetPosition().y();
179
180 fInitGamma = x * px + y * py;
181
182 fInitGamma =
183 fInitGamma / std::sqrt(px * px + py * py) / std::sqrt(x * x + y * y);
184
185 fInitGamma = std::acos(fInitGamma * rad);
186
187 if(fInitGamma / deg > 90.0)
188 {
189 fInitGamma = 180 * deg - fInitGamma;
190 }
191 }
192 // Record Photons that missed the photon detector but escaped from readout
193 if(!thePostPV && trackInformation->IsStatus(EscapedFromReadOut))
194 {
195 //G4cout << "SteppingAction: status = EscapedFromReadOut" << G4endl;
197 // UpdateHistogramSuccess(thePostPoint,theTrack);
199
200 return;
201 }
202
203 // Assumed photons are originated at the fiber OR
204 // the fiber is the first material the photon hits
205 switch(theStatus)
206 {
207 // Exiting the fiber
208 case FresnelRefraction:
209 case SameMaterial:
211
212 if(thePostPVname == "WLSFiber" || thePostPVname == "Clad1" ||
213 thePostPVname == "Clad2")
214 {
215 if(trackInformation->IsStatus(OutsideOfFiber))
216 trackInformation->AddStatusFlag(InsideOfFiber);
217
218 // Set the Exit flag when the photon refracted out of the fiber
219 }
220 else if(trackInformation->IsStatus(InsideOfFiber))
221 {
222 // EscapedFromReadOut if the z position is the same as fiber's end
223 if(theTrack->GetPosition().z() == fDetector->GetWLSFiberEnd())
224 {
225 trackInformation->AddStatusFlag(EscapedFromReadOut);
226 fCounterEnd++;
228 }
229 else // Escaped from side
230 {
231 trackInformation->AddStatusFlag(EscapedFromSide);
232 trackInformation->SetExitPosition(theTrack->GetPosition());
233 // UpdateHistogramEscape(thePostPoint,theTrack);
234
235 fCounterMid++;
238 }
239
240 trackInformation->AddStatusFlag(OutsideOfFiber);
241 trackInformation->SetExitPosition(theTrack->GetPosition());
242 }
243
244 return;
245
246 // Internal Reflections
247 case TotalInternalReflection:
248
250
251 // Kill the track if it's number of bounces exceeded the limit
253 {
254 theTrack->SetTrackStatus(fStopAndKill);
255 trackInformation->AddStatusFlag(murderee);
257 G4cout << "\n Bounce Limit Exceeded" << G4endl;
258 return;
259 }
260 break;
261
262 case FresnelReflection:
263
266
267 if(thePrePVname == "WLSFiber")
268 {
271 }
272 else if(thePrePVname == "Clad1")
273 {
276 }
277 else if(thePrePVname == "Clad2")
278 {
281 }
282
283 // Determine if the photon has reflected off the read-out end
284 if(theTrack->GetPosition().z() == fDetector->GetWLSFiberEnd())
285 {
286 if(!trackInformation->IsStatus(ReflectedAtReadOut) &&
287 trackInformation->IsStatus(InsideOfFiber))
288 {
289 trackInformation->AddStatusFlag(ReflectedAtReadOut);
290
292 theStatus == TotalInternalReflection)
293 {
294 theTrack->SetTrackStatus(fStopAndKill);
295 trackInformation->AddStatusFlag(murderee);
296 // UpdateHistogramReflect(thePostPoint,theTrack);
298 return;
299 }
300 }
301 }
302 return;
303
304 // Reflection off the mirror
305 case LambertianReflection:
306 case LobeReflection:
307 case SpikeReflection:
308
310 // Check if it hits the mirror
311 if(thePostPVname == "Mirror")
312 {
313 trackInformation->AddStatusFlag(ReflectedAtMirror);
315 }
316 return;
317
318 // Detected by a detector
319 case Detection:
320 // Detected automatically with G4OpBoundaryProcess->InvokeSD set true
321
322 // Stop Tracking when it hits the detector's surface
324 theTrack->SetTrackStatus(fStopAndKill);
325 return;
326
327 default:
328 break;
329 }
330
331 // Check for absorbed photons
332 if(theTrack->GetTrackStatus() != fAlive &&
333 trackInformation->IsStatus(InsideOfFiber))
334 {
335 // UpdateHistogramAbsorb(thePostPoint,theTrack);
337 return;
338 }
339}
Definition of the WLSDetectorConstruction class.
Definition of the WLSEventAction class.
Definition of the WLSPhotonDetSD class.
Definition of the WLSSteppingActionMessenger class.
Definition of the WLSSteppingAction class.
Definition of the WLSUserTrackInformation class.
WLSDetectorConstruction * fDetector
G4OpBoundaryProcess * fOpProcess
WLSSteppingActionMessenger * fSteppingMessenger
void UserSteppingAction(const G4Step *) override
~WLSSteppingAction() override
WLSEventAction * fEventAction
WLSSteppingAction(WLSDetectorConstruction *, WLSEventAction *)

Applications | User Support | Publications | Collaboration