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 // class G4RegularNavigation 27 // 28 // Class description: 29 // 30 // Utility for fast navigation in volumes containing a regular 31 // parameterisation. If two contiguous voxels have the same material, 32 // navigation does not stop at the surface 33 34 // History: 35 // - Created. P. Arce, May 2007 36 // -------------------------------------------------------------------- 37 #ifndef G4RegularNavigation_HH 38 #define G4RegularNavigation_HH 39 40 #include <vector> 41 42 #include "G4Types.hh" 43 #include "G4ThreeVector.hh" 44 #include "G4VNavigation.hh" 45 46 class G4NormalNavigation; 47 class G4VPhysicalVolume; 48 class G4Navigator; 49 class G4NavigationHistory; 50 51 class G4RegularNavigation : public G4VNavigation 52 { 53 public: // with description 54 55 G4RegularNavigation(); 56 ~G4RegularNavigation() override; 57 58 G4bool LevelLocate( G4NavigationHistory& history, 59 const G4VPhysicalVolume* blockedVol, 60 const G4int blockedNum, 61 const G4ThreeVector& globalPoint, 62 const G4ThreeVector* globalDirection, 63 const G4bool pLocatedOnEdge, 64 G4ThreeVector& localPoint ) final; 65 // Locate point using its position with respect to regular 66 // parameterisation container volume. 67 68 G4double ComputeStep( const G4ThreeVector& globalPoint, 69 const G4ThreeVector& globalDirection, 70 const G4double currentProposedStepLength, 71 G4double& newSafety, 72 G4NavigationHistory& history, 73 G4bool& validExitNormal, 74 G4ThreeVector& exitNormal, 75 G4bool& exiting, 76 G4bool& entering, 77 G4VPhysicalVolume *(*pBlockedPhysical), 78 G4int& blockedReplicaNo ) final; 79 // Method never called because to be called the daughter has to be a 80 // 'regular' volume. This would only happen if the track is in the 81 // mother of voxels volume. But the voxels fill completely their mother, 82 // so when a track enters the mother it automatically enters a voxel. 83 84 G4double ComputeStepSkippingEqualMaterials( 85 G4ThreeVector& localPoint, 86 const G4ThreeVector& globalDirection, 87 const G4double currentProposedStepLength, 88 G4double& newSafety, 89 G4NavigationHistory& history, 90 G4bool& validExitNormal, 91 G4ThreeVector& exitNormal, 92 G4bool& exiting, 93 G4bool& entering, 94 G4VPhysicalVolume *(*pBlockedPhysical), 95 G4int& blockedReplicaNo, 96 G4VPhysicalVolume* pCurrentPhysical); 97 // Compute the step skipping surfaces when they separate voxels with 98 // equal materials. Loop to voxels until a different material is found: 99 // invokes G4NormalNavigation::ComputeStep() in each voxel and move the 100 // point to the next voxel. 101 102 G4double ComputeSafety( const G4ThreeVector& localPoint, 103 const G4NavigationHistory& history, 104 const G4double pProposedMaxLength = DBL_MAX ) final; 105 // Method never called because to be called the daughter has to be a 106 // 'regular' volume. This would only happen if the track is in the 107 // mother of voxels volume. But the voxels fill completely their mother, 108 // so when a track enters the mother it automatically enters a voxel. 109 110 public: // without description 111 112 void SetNormalNavigation( G4NormalNavigation* fnormnav ) 113 { fnormalNav = fnormnav; } 114 115 private: 116 117 G4NormalNavigation* fnormalNav = nullptr; 118 G4double kCarTolerance; 119 G4double fMinStep; 120 121 G4bool fLastStepWasZero = false; 122 // Whether the last ComputeStep moved Zero. Used to check for edges. 123 G4int fNumberZeroSteps = 0; 124 // Number of preceding moves that were Zero. Reset to 0 after finite step 125 G4int fActionThreshold_NoZeroSteps = 2; 126 // After this many failed/zero steps, act (push etc) 127 G4int fAbandonThreshold_NoZeroSteps = 25; 128 // After this many failed/zero steps, abandon track 129 G4int fNoStepsAllowed = 10000; 130 // Maximum number of steps a track can travel skipping voxels (if there are more, track is assumed to be stuck and it is killed) 131 }; 132 133 #endif 134