Loading...
Searching...
No Matches
Public Member Functions | Private Member Functions | Private Attributes | List of all members
WLSSteppingAction Class Reference

#include <Doxymodules_optical.h>

Inheritance diagram for WLSSteppingAction:
G4UserSteppingAction

Public Member Functions

 WLSSteppingAction (WLSDetectorConstruction *, WLSEventAction *)
 
 ~WLSSteppingAction () override
 
void UserSteppingAction (const G4Step *) override
 
void SetBounceLimit (G4int)
 
G4int GetNumberOfBounces ()
 
G4int GetNumberOfClad1Bounces ()
 
G4int GetNumberOfClad2Bounces ()
 
G4int GetNumberOfWLSBounces ()
 
G4int ResetSuccessCounter ()
 

Private Member Functions

void ResetCounters ()
 

Private Attributes

G4int fBounceLimit = 100000
 
G4int fCounterEnd = 0
 
G4int fCounterMid = 0
 
G4int fCounterBounce = 0
 
G4int fCounterWLSBounce = 0
 
G4int fCounterClad1Bounce = 0
 
G4int fCounterClad2Bounce = 0
 
G4double fInitGamma = -1.
 
G4double fInitTheta = -1.
 
G4OpBoundaryProcessfOpProcess = nullptr
 
WLSDetectorConstructionfDetector = nullptr
 
WLSSteppingActionMessengerfSteppingMessenger = nullptr
 
WLSEventActionfEventAction = nullptr
 

Detailed Description

Definition at line 102 of file Doxymodules_optical.h.

Constructor & Destructor Documentation

◆ WLSSteppingAction()

WLSSteppingAction::WLSSteppingAction ( WLSDetectorConstruction detector,
WLSEventAction event 
)

Definition at line 58 of file WLSSteppingAction.cc.

60 : fDetector(detector)
61 , fEventAction(event)
62{
65}
WLSDetectorConstruction * fDetector
WLSSteppingActionMessenger * fSteppingMessenger
WLSEventAction * fEventAction

◆ ~WLSSteppingAction()

WLSSteppingAction::~WLSSteppingAction ( )
override

Definition at line 69 of file WLSSteppingAction.cc.

69{ delete fSteppingMessenger; }

Member Function Documentation

◆ UserSteppingAction()

void WLSSteppingAction::UserSteppingAction ( const G4Step theStep)
override

Definition at line 108 of file WLSSteppingAction.cc.

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}
G4OpBoundaryProcess * fOpProcess

◆ SetBounceLimit()

void WLSSteppingAction::SetBounceLimit ( G4int  i)

Definition at line 73 of file WLSSteppingAction.cc.

73{ fBounceLimit = i; }

◆ GetNumberOfBounces()

G4int WLSSteppingAction::GetNumberOfBounces ( )

Definition at line 77 of file WLSSteppingAction.cc.

77{ return fCounterBounce; }

◆ GetNumberOfClad1Bounces()

G4int WLSSteppingAction::GetNumberOfClad1Bounces ( )

Definition at line 81 of file WLSSteppingAction.cc.

82{
84}

◆ GetNumberOfClad2Bounces()

G4int WLSSteppingAction::GetNumberOfClad2Bounces ( )

Definition at line 88 of file WLSSteppingAction.cc.

89{
91}

◆ GetNumberOfWLSBounces()

G4int WLSSteppingAction::GetNumberOfWLSBounces ( )

Definition at line 95 of file WLSSteppingAction.cc.

95{ return fCounterWLSBounce; }

◆ ResetSuccessCounter()

G4int WLSSteppingAction::ResetSuccessCounter ( )

Definition at line 99 of file WLSSteppingAction.cc.

100{
101 G4int temp = fCounterEnd;
102 fCounterEnd = 0;
103 return temp;
104}

◆ ResetCounters()

void WLSSteppingAction::ResetCounters ( )
inlineprivate

Definition at line 93 of file WLSSteppingAction.hh.

94 {
99 fInitGamma = -1;
100 fInitTheta = -1;
101 }

Member Data Documentation

◆ fBounceLimit

G4int WLSSteppingAction::fBounceLimit = 100000
private

Definition at line 68 of file WLSSteppingAction.hh.

◆ fCounterEnd

G4int WLSSteppingAction::fCounterEnd = 0
private

Definition at line 70 of file WLSSteppingAction.hh.

◆ fCounterMid

G4int WLSSteppingAction::fCounterMid = 0
private

Definition at line 72 of file WLSSteppingAction.hh.

◆ fCounterBounce

G4int WLSSteppingAction::fCounterBounce = 0
private

Definition at line 74 of file WLSSteppingAction.hh.

◆ fCounterWLSBounce

G4int WLSSteppingAction::fCounterWLSBounce = 0
private

Definition at line 76 of file WLSSteppingAction.hh.

◆ fCounterClad1Bounce

G4int WLSSteppingAction::fCounterClad1Bounce = 0
private

Definition at line 78 of file WLSSteppingAction.hh.

◆ fCounterClad2Bounce

G4int WLSSteppingAction::fCounterClad2Bounce = 0
private

Definition at line 80 of file WLSSteppingAction.hh.

◆ fInitGamma

G4double WLSSteppingAction::fInitGamma = -1.
private

Definition at line 83 of file WLSSteppingAction.hh.

◆ fInitTheta

G4double WLSSteppingAction::fInitTheta = -1.
private

Definition at line 85 of file WLSSteppingAction.hh.

◆ fOpProcess

G4OpBoundaryProcess* WLSSteppingAction::fOpProcess = nullptr
private

Definition at line 87 of file WLSSteppingAction.hh.

◆ fDetector

WLSDetectorConstruction* WLSSteppingAction::fDetector = nullptr
private

Definition at line 89 of file WLSSteppingAction.hh.

◆ fSteppingMessenger

WLSSteppingActionMessenger* WLSSteppingAction::fSteppingMessenger = nullptr
private

Definition at line 90 of file WLSSteppingAction.hh.

◆ fEventAction

WLSEventAction* WLSSteppingAction::fEventAction = nullptr
private

Definition at line 91 of file WLSSteppingAction.hh.


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

Applications | User Support | Publications | Collaboration