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 // G4SteppingManager 27 // 28 // Class description: 29 // 30 // This is the class which plays an essential role in tracking particles. 31 // It takes cares of all message passing between objects in the different 32 // categories (for example, geometry - including transportation, interactions 33 // in matter, etc). Its public method 'stepping' steers to step the particle. 34 // Only used within the Geant4 kernel. 35 36 // Contact: 37 // Questions and comments to this code should be sent to 38 // Katsuya Amako (e-mail: Katsuya.Amako@kek.jp) 39 // Takashi Sasaki (e-mail: Takashi.Sasaki@kek.jp) 40 // 41 // History: 42 // 12.03.1998, H.Kurashige - modified for use of G4ParticleChange 43 // -------------------------------------------------------------------- 44 #ifndef G4SteppingManager_hh 45 #define G4SteppingManager_hh 1 46 47 #include "G4LogicalVolume.hh" // Include from 'geometry' 48 #include "G4Navigator.hh" // Include from 'geometry' 49 #include "G4NoProcess.hh" // Include from 'processes' 50 #include "G4ProcessManager.hh" // Include from 'processes' 51 #include "G4Step.hh" // Include from 'tracking' 52 #include "G4StepPoint.hh" // Include from 'tracking' 53 #include "G4StepStatus.hh" // Include from 'tracking' 54 #include "G4TouchableHandle.hh" // Include from 'geometry' 55 #include "G4Track.hh" // Include from 'tracking' 56 #include "G4TrackStatus.hh" // Include from 'tracking' 57 #include "G4TrackVector.hh" // Include from 'tracking' 58 #include "G4UserSteppingAction.hh" // Include from 'tracking' 59 #include "G4VPhysicalVolume.hh" // Include from 'geometry' 60 #include "G4VSteppingVerbose.hh" // Include from 'tracking' 61 #include "Randomize.hh" // Include from 'global' 62 #include "globals.hh" // Include from 'global' 63 64 #include <iomanip> // Include from 'system' 65 #include <vector> // Include from 'system' 66 67 using G4SelectedAtRestDoItVector = std::vector<G4int>; 68 using G4SelectedAlongStepDoItVector = std::vector<G4int>; 69 using G4SelectedPostStepDoItVector = std::vector<G4int>; 70 71 class G4VSensitiveDetector; 72 73 class G4SteppingManager 74 { 75 public: 76 // Constructor/Destructor 77 78 // SteppingManger should be dynamically allocated, therefore 79 // you need to invoke new() when you call this constructor. 80 // "Secodary track vector" will be dynamically created by this 81 // constructor. G4UserSteppingAction will be also created 82 // in this constructor, and "this" pointer will be passed to 83 // G4UserSteppingAction. 84 G4SteppingManager(); 85 ~G4SteppingManager(); 86 87 // Get/Set functions 88 89 const G4TrackVector* GetSecondary() const; 90 void SetUserAction(G4UserSteppingAction* apAction); 91 G4Track* GetTrack() const; 92 void SetVerboseLevel(G4int vLevel); 93 void SetVerbose(G4VSteppingVerbose*); 94 G4Step* GetStep() const; 95 void SetNavigator(G4Navigator* value); 96 97 // Other member functions 98 99 // Steers to move the give particle from the TrackingManger by one Step. 100 G4StepStatus Stepping(); 101 102 // Sets up initial track information (enegry, position, etc) to 103 // the PreStepPoint of the G4Step. This funciton has to be called 104 // just once before the stepping loop in the "TrackingManager". 105 void SetInitialStep(G4Track* valueTrack); 106 107 // Get methods 108 void GetProcessNumber(); 109 G4double GetPhysicalStep(); 110 G4double GetGeometricalStep(); 111 G4double GetCorrectedStep(); 112 G4bool GetPreStepPointIsGeom(); 113 G4bool GetFirstStep(); 114 G4StepStatus GetfStepStatus(); 115 G4double GetTempInitVelocity(); 116 G4double GetTempVelocity(); 117 G4double GetMass(); 118 G4double GetsumEnergyChange(); 119 G4VParticleChange* GetfParticleChange(); 120 G4Track* GetfTrack(); 121 G4TrackVector* GetfSecondary(); 122 G4Step* GetfStep(); 123 G4StepPoint* GetfPreStepPoint(); 124 G4StepPoint* GetfPostStepPoint(); 125 G4VPhysicalVolume* GetfCurrentVolume(); 126 G4VSensitiveDetector* GetfSensitive(); 127 G4VProcess* GetfCurrentProcess(); 128 G4ProcessVector* GetfAtRestDoItVector(); 129 G4ProcessVector* GetfAlongStepDoItVector(); 130 G4ProcessVector* GetfPostStepDoItVector(); 131 G4ProcessVector* GetfAlongStepGetPhysIntVector(); 132 G4ProcessVector* GetfPostStepGetPhysIntVector(); 133 G4ProcessVector* GetfAtRestGetPhysIntVector(); 134 G4double GetcurrentMinimumStep(); 135 G4double GetnumberOfInteractionLengthLeft(); 136 std::size_t GetfAtRestDoItProcTriggered(); 137 std::size_t GetfAlongStepDoItProcTriggered(); 138 std::size_t GetfPostStepDoItProcTriggered(); 139 G4int GetfN2ndariesAtRestDoIt(); 140 G4int GetfN2ndariesAlongStepDoIt(); 141 G4int GetfN2ndariesPostStepDoIt(); 142 G4Navigator* GetfNavigator(); 143 G4int GetverboseLevel(); 144 std::size_t GetMAXofAtRestLoops(); 145 std::size_t GetMAXofAlongStepLoops(); 146 std::size_t GetMAXofPostStepLoops(); 147 G4SelectedAtRestDoItVector* GetfSelectedAtRestDoItVector(); 148 G4SelectedAlongStepDoItVector* GetfSelectedAlongStepDoItVector(); 149 G4SelectedPostStepDoItVector* GetfSelectedPostStepDoItVector(); 150 G4double GetfPreviousStepSize(); 151 const G4TouchableHandle& GetTouchableHandle(); 152 G4SteppingControl GetStepControlFlag(); 153 G4UserSteppingAction* GetUserAction(); 154 G4double GetphysIntLength(); 155 G4ForceCondition GetfCondition(); 156 G4GPILSelection GetfGPILSelection(); 157 158 private: 159 // Member functions 160 161 // Calculate corresponding physical length from the mean free path 162 // left for each discrete physics process. The minimum allowable 163 // step for each continuous process will be also calculated. 164 void DefinePhysicalStepLength(); 165 166 void InvokeAtRestDoItProcs(); 167 void InvokeAlongStepDoItProcs(); 168 void InvokePostStepDoItProcs(); 169 void InvokePSDIP(size_t); // 170 G4int ProcessSecondariesFromParticleChange(); 171 172 // Return the estimated safety value at the PostStepPoint 173 G4double CalculateSafety(); 174 175 // Member data 176 177 static const size_t SizeOfSelectedDoItVector = 100; 178 179 G4bool KillVerbose = false; 180 181 G4UserSteppingAction* fUserSteppingAction = nullptr; 182 183 G4VSteppingVerbose* fVerbose = nullptr; 184 185 G4double PhysicalStep = 0.0; 186 G4double GeometricalStep = 0.0; 187 G4double CorrectedStep = 0.0; 188 G4bool PreStepPointIsGeom = false; 189 G4bool FirstStep = false; 190 G4StepStatus fStepStatus = fUndefined; 191 192 G4double TempInitVelocity = 0.0; 193 G4double TempVelocity = 0.0; 194 G4double Mass = 0.0; 195 196 G4double sumEnergyChange = 0.0; 197 198 G4VParticleChange* fParticleChange = nullptr; 199 G4Track* fTrack = nullptr; 200 G4TrackVector* fSecondary = nullptr; 201 G4Step* fStep = nullptr; 202 G4StepPoint* fPreStepPoint = nullptr; 203 G4StepPoint* fPostStepPoint = nullptr; 204 205 G4VPhysicalVolume* fCurrentVolume = nullptr; 206 G4VSensitiveDetector* fSensitive = nullptr; 207 G4VProcess* fCurrentProcess = nullptr; // Pointer to process of which DoIt() or 208 // GetPhysicalInteractionLength() has been just executed. 209 210 G4ProcessVector* fAtRestDoItVector = nullptr; 211 G4ProcessVector* fAlongStepDoItVector = nullptr; 212 G4ProcessVector* fPostStepDoItVector = nullptr; 213 214 G4ProcessVector* fAtRestGetPhysIntVector = nullptr; 215 G4ProcessVector* fAlongStepGetPhysIntVector = nullptr; 216 G4ProcessVector* fPostStepGetPhysIntVector = nullptr; 217 218 std::size_t MAXofAtRestLoops = 0; 219 std::size_t MAXofAlongStepLoops = 0; 220 std::size_t MAXofPostStepLoops = 0; 221 222 std::size_t fAtRestDoItProcTriggered = 0; 223 std::size_t fAlongStepDoItProcTriggered = 0; 224 std::size_t fPostStepDoItProcTriggered = 0; 225 226 G4int fN2ndariesAtRestDoIt = 0; 227 G4int fN2ndariesAlongStepDoIt = 0; 228 G4int fN2ndariesPostStepDoIt = 0; 229 // These are the numbers of secondaries generated by the process 230 // just executed. 231 232 G4Navigator* fNavigator = nullptr; 233 234 G4int verboseLevel = 0; 235 236 G4SelectedAtRestDoItVector* fSelectedAtRestDoItVector = nullptr; 237 G4SelectedAlongStepDoItVector* fSelectedAlongStepDoItVector = nullptr; 238 G4SelectedPostStepDoItVector* fSelectedPostStepDoItVector = nullptr; 239 240 G4double fPreviousStepSize = 0.0; 241 242 G4TouchableHandle fTouchableHandle; 243 244 G4SteppingControl StepControlFlag = NormalCondition; 245 246 G4double kCarTolerance = 0.0; // Cached geometrical tolerance on surface 247 G4double proposedSafety = 0.0; // This keeps the minimum safety value proposed by AlongStepGPILs. 248 G4ThreeVector endpointSafOrigin; 249 G4double endpointSafety = 0.0; // To get the true safety value at the PostStepPoint, you have to 250 // subtract the distance to 'endpointSafOrigin' from this value. 251 G4double physIntLength = 0.0; 252 G4ForceCondition fCondition = InActivated; 253 G4GPILSelection fGPILSelection = NotCandidateForSelection; 254 // Above three variables are for the method 255 // DefinePhysicalStepLength(). To pass these information to 256 // the method Verbose, they are kept at here. Need a more 257 // elegant mechanism. 258 259 G4NoProcess const* fNoProcess = nullptr; // Used in InvokeAtRestDoItProcs() method to flag the 260 // process of any stable ion at rest. 261 }; 262 263 //******************************************************************* 264 // 265 // Inline functions 266 // 267 //******************************************************************* 268 269 inline G4double G4SteppingManager::GetPhysicalStep() { return PhysicalStep; } 270 271 inline G4double G4SteppingManager::GetGeometricalStep() { return GeometricalStep; } 272 273 inline G4double G4SteppingManager::GetCorrectedStep() { return CorrectedStep; } 274 275 inline G4bool G4SteppingManager::GetPreStepPointIsGeom() { return PreStepPointIsGeom; } 276 277 inline G4bool G4SteppingManager::GetFirstStep() { return FirstStep; } 278 279 inline G4StepStatus G4SteppingManager::GetfStepStatus() { return fStepStatus; } 280 281 inline G4double G4SteppingManager::GetTempInitVelocity() { return TempInitVelocity; } 282 283 inline G4double G4SteppingManager::GetTempVelocity() { return TempVelocity; } 284 285 inline G4double G4SteppingManager::GetMass() { return Mass; } 286 287 inline G4double G4SteppingManager::GetsumEnergyChange() { return sumEnergyChange; } 288 289 inline G4VParticleChange* G4SteppingManager::GetfParticleChange() { return fParticleChange; } 290 291 inline G4Track* G4SteppingManager::GetfTrack() { return fTrack; } 292 293 inline G4TrackVector* G4SteppingManager::GetfSecondary() { return fStep->GetfSecondary(); } 294 295 inline G4Step* G4SteppingManager::GetfStep() { return fStep; } 296 297 inline G4StepPoint* G4SteppingManager::GetfPreStepPoint() { return fPreStepPoint; } 298 299 inline G4StepPoint* G4SteppingManager::GetfPostStepPoint() { return fPostStepPoint; } 300 301 inline G4VPhysicalVolume* G4SteppingManager::GetfCurrentVolume() { return fCurrentVolume; } 302 303 inline G4VSensitiveDetector* G4SteppingManager::GetfSensitive() { return fSensitive; } 304 305 inline G4VProcess* G4SteppingManager::GetfCurrentProcess() { return fCurrentProcess; } 306 307 inline G4ProcessVector* G4SteppingManager::GetfAtRestDoItVector() { return fAtRestDoItVector; } 308 309 inline G4ProcessVector* G4SteppingManager::GetfAlongStepDoItVector() 310 { 311 return fAlongStepDoItVector; 312 } 313 314 inline G4ProcessVector* G4SteppingManager::GetfPostStepDoItVector() { return fPostStepDoItVector; } 315 316 inline G4ProcessVector* G4SteppingManager::GetfAtRestGetPhysIntVector() 317 { 318 return fAtRestGetPhysIntVector; 319 } 320 321 inline G4ProcessVector* G4SteppingManager::GetfAlongStepGetPhysIntVector() 322 { 323 return fAlongStepGetPhysIntVector; 324 } 325 326 inline G4ProcessVector* G4SteppingManager::GetfPostStepGetPhysIntVector() 327 { 328 return fPostStepGetPhysIntVector; 329 } 330 331 inline size_t G4SteppingManager::GetMAXofAtRestLoops() { return MAXofAtRestLoops; } 332 333 inline size_t G4SteppingManager::GetMAXofAlongStepLoops() { return MAXofAlongStepLoops; } 334 335 inline size_t G4SteppingManager::GetMAXofPostStepLoops() { return MAXofPostStepLoops; } 336 337 inline size_t G4SteppingManager::GetfAtRestDoItProcTriggered() { return fAtRestDoItProcTriggered; } 338 339 inline size_t G4SteppingManager::GetfAlongStepDoItProcTriggered() 340 { 341 return fAtRestDoItProcTriggered; 342 } 343 344 inline size_t G4SteppingManager::GetfPostStepDoItProcTriggered() 345 { 346 return fPostStepDoItProcTriggered; 347 } 348 349 inline G4int G4SteppingManager::GetfN2ndariesAtRestDoIt() { return fN2ndariesAtRestDoIt; } 350 351 inline G4int G4SteppingManager::GetfN2ndariesAlongStepDoIt() { return fN2ndariesAlongStepDoIt; } 352 353 inline G4int G4SteppingManager::GetfN2ndariesPostStepDoIt() { return fN2ndariesPostStepDoIt; } 354 355 inline G4Navigator* G4SteppingManager::GetfNavigator() { return fNavigator; } 356 357 inline G4int G4SteppingManager::GetverboseLevel() { return verboseLevel; } 358 359 inline G4SelectedAtRestDoItVector* G4SteppingManager::GetfSelectedAtRestDoItVector() 360 { 361 return fSelectedAtRestDoItVector; 362 } 363 364 inline G4SelectedAlongStepDoItVector* G4SteppingManager::GetfSelectedAlongStepDoItVector() 365 { 366 return fSelectedAlongStepDoItVector; 367 } 368 369 inline G4SelectedPostStepDoItVector* G4SteppingManager::GetfSelectedPostStepDoItVector() 370 { 371 return fSelectedPostStepDoItVector; 372 } 373 374 inline G4double G4SteppingManager::GetfPreviousStepSize() { return fPreviousStepSize; } 375 376 inline const G4TouchableHandle& G4SteppingManager::GetTouchableHandle() { return fTouchableHandle; } 377 378 inline G4SteppingControl G4SteppingManager::GetStepControlFlag() { return StepControlFlag; } 379 380 inline G4double G4SteppingManager::GetphysIntLength() { return physIntLength; } 381 382 inline G4ForceCondition G4SteppingManager::GetfCondition() { return fCondition; } 383 384 inline G4GPILSelection G4SteppingManager::GetfGPILSelection() { return fGPILSelection; } 385 386 inline const G4TrackVector* G4SteppingManager::GetSecondary() const 387 { 388 return fStep->GetSecondary(); 389 } 390 391 inline void G4SteppingManager::SetNavigator(G4Navigator* value) { fNavigator = value; } 392 393 inline void G4SteppingManager::SetUserAction(G4UserSteppingAction* apAction) 394 { 395 fUserSteppingAction = apAction; 396 } 397 398 inline G4UserSteppingAction* G4SteppingManager::GetUserAction() { return fUserSteppingAction; } 399 400 inline G4Track* G4SteppingManager::GetTrack() const { return fTrack; } 401 402 inline void G4SteppingManager::SetVerboseLevel(G4int vLevel) { verboseLevel = vLevel; } 403 404 inline void G4SteppingManager::SetVerbose(G4VSteppingVerbose* yourVerbose) 405 { 406 fVerbose = yourVerbose; 407 } 408 409 inline G4Step* G4SteppingManager::GetStep() const { return fStep; } 410 411 inline G4double G4SteppingManager::CalculateSafety() 412 { 413 return std::max( 414 endpointSafety - (endpointSafOrigin - fPostStepPoint->GetPosition()).mag(), kCarTolerance); 415 } 416 417 #endif 418