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