|
||||
File indexing completed on 2025-01-18 09:59:04
0001 // 0002 // ******************************************************************** 0003 // * License and Disclaimer * 0004 // * * 0005 // * The Geant4 software is copyright of the Copyright Holders of * 0006 // * the Geant4 Collaboration. It is provided under the terms and * 0007 // * conditions of the Geant4 Software License, included in the file * 0008 // * LICENSE and available at http://cern.ch/geant4/license . These * 0009 // * include a list of copyright holders. * 0010 // * * 0011 // * Neither the authors of this software system, nor their employing * 0012 // * institutes,nor the agencies providing financial support for this * 0013 // * work make any representation or warranty, express or implied, * 0014 // * regarding this software system or assume any liability for its * 0015 // * use. Please see the license in the file LICENSE and URL above * 0016 // * for the full disclaimer and the limitation of liability. * 0017 // * * 0018 // * This code implementation is the result of the scientific and * 0019 // * technical work of the GEANT4 collaboration. * 0020 // * By using, copying, modifying or distributing the software (or * 0021 // * any work based on the software) you agree to acknowledge its * 0022 // * use in resulting scientific publications, and indicate your * 0023 // * acceptance of all terms of the Geant4 Software license. * 0024 // ******************************************************************** 0025 // 0026 // class G4SafetyHelper 0027 // 0028 // Class description: 0029 // 0030 // This class is a helper for physics processes which require 0031 // knowledge of the safety, and the step size for the 'mass' geometry 0032 0033 // First version: J.Apostolakis, July 5th, 2006 0034 // -------------------------------------------------------------------- 0035 #ifndef G4SAFETYHELPER_HH 0036 #define G4SAFETYHELPER_HH 1 0037 0038 #include <vector> 0039 0040 #include "G4Types.hh" 0041 #include "G4ThreeVector.hh" 0042 #include "G4Navigator.hh" 0043 0044 class G4PathFinder; 0045 0046 class G4SafetyHelper 0047 { 0048 public: // with description 0049 0050 G4SafetyHelper(); 0051 ~G4SafetyHelper(); 0052 // Constructor and destructor 0053 0054 G4double CheckNextStep( const G4ThreeVector& position, 0055 const G4ThreeVector& direction, 0056 const G4double currentMaxStep, 0057 G4double& newSafety ); 0058 // Return linear step for mass geometry 0059 0060 G4double ComputeSafety( const G4ThreeVector& pGlobalPoint, 0061 G4double maxRadius = DBL_MAX ); 0062 // Return safety for all geometries. 0063 // 0064 // The 2nd argument is the radius of your interest (e.g. maximum 0065 // displacement). Giving this you can reduce the average computational 0066 // cost. If the second argument is not given, this is the real 0067 // isotropic safety 0068 0069 void Locate(const G4ThreeVector& pGlobalPoint, 0070 const G4ThreeVector& direction); 0071 // Locate the point for all geometries 0072 0073 void ReLocateWithinVolume(const G4ThreeVector& pGlobalPoint ); 0074 // Relocate the point in the volume of interest 0075 0076 inline void EnableParallelNavigation(G4bool parallel); 0077 // To have parallel worlds considered, must be true. 0078 // Alternative is to use single (mass) Navigator directly 0079 0080 void InitialiseNavigator(); 0081 // Check for new navigator for tracking, and reinitialise pointer 0082 0083 inline G4int SetVerboseLevel( G4int lev ); 0084 inline G4VPhysicalVolume* GetWorldVolume(); 0085 inline void SetCurrentSafety(G4double val, const G4ThreeVector& pos); 0086 0087 public: // without description 0088 0089 void InitialiseHelper(); 0090 0091 private: 0092 0093 G4PathFinder* fpPathFinder = nullptr; 0094 G4Navigator* fpMassNavigator = nullptr; 0095 0096 G4bool fUseParallelGeometries = false; 0097 // Flag whether to use PathFinder or single (mass) Navigator directly 0098 // By default, one geometry only 0099 G4bool fFirstCall = true; 0100 // Flag of first call 0101 G4int fVerbose = 0; 0102 // Whether to print warning in case of move outside safety 0103 0104 // State used during tracking -- for optimisation 0105 0106 G4ThreeVector fLastSafetyPosition; 0107 G4double fLastSafety = 0.0; 0108 0109 // const G4double fRecomputeFactor = 0.0; 0110 // parameter for further optimisation: 0111 // if ( move < fact*safety ) do fast recomputation of safety 0112 0113 // End State (tracking) 0114 }; 0115 0116 // -------------------------------------------------------------------- 0117 // Inline definitions 0118 // -------------------------------------------------------------------- 0119 0120 inline G4int G4SafetyHelper::SetVerboseLevel( G4int lev ) 0121 { 0122 G4int oldlv = fVerbose; 0123 fVerbose = lev; 0124 return oldlv; 0125 } 0126 0127 inline 0128 void G4SafetyHelper::EnableParallelNavigation(G4bool parallel) 0129 { 0130 fUseParallelGeometries = parallel; 0131 } 0132 0133 inline 0134 G4VPhysicalVolume* G4SafetyHelper::GetWorldVolume() 0135 { 0136 return fpMassNavigator->GetWorldVolume(); 0137 } 0138 0139 inline 0140 void G4SafetyHelper::SetCurrentSafety(G4double val, const G4ThreeVector& pos) 0141 { 0142 fLastSafety = val; 0143 fLastSafetyPosition = pos; 0144 } 0145 0146 #endif
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |