File indexing completed on 2025-01-18 09:57:53
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033
0034
0035
0036
0037
0038
0039
0040
0041
0042
0043
0044
0045
0046 #ifndef G4AllITFinder_h
0047 #define G4AllITFinder_h 1
0048 #include "globals.hh"
0049 #include "G4ITType.hh"
0050 #include "G4ThreeVector.hh"
0051 #include "G4KDTreeResult.hh"
0052
0053 #include <map>
0054 #include <vector>
0055 #include <memory>
0056
0057 class G4IT;
0058 class G4VITFinder;
0059 class G4ITBox;
0060 class G4Track;
0061 template<typename T>
0062 class G4ITFinder;
0063
0064
0065
0066
0067
0068
0069
0070 class G4AllITFinder
0071 {
0072 public:
0073 static G4AllITFinder* Instance();
0074 static void DeleteInstance();
0075
0076
0077
0078
0079
0080 ~G4AllITFinder();
0081
0082 G4VITFinder* GetInstance(G4ITType);
0083 G4ITBox* GetBox(const G4Track*);
0084
0085 void RegisterManager(G4VITFinder* manager);
0086 void Push(G4Track* track);
0087
0088
0089
0090
0091
0092 void SetVerboseLevel(G4int level)
0093 {
0094 fVerbose = level;
0095 }
0096 G4int GetVerboseLevel()
0097 {
0098 return fVerbose;
0099 }
0100
0101 void UpdatePositionMap();
0102
0103 template<typename T>
0104 inline G4KDTreeResultHandle FindNearest(const G4ThreeVector& pos,
0105 const T* it);
0106 template<typename T>
0107 inline G4KDTreeResultHandle FindNearest(const T* it0, const T* it);
0108 template<typename T>
0109 inline G4KDTreeResultHandle FindNearestInRange(const G4ThreeVector& pos,
0110 const T* it,
0111 G4double range);
0112 template<typename T>
0113 inline G4KDTreeResultHandle FindNearestInRange(const T* it0,
0114 const T* it,
0115 G4double range);
0116
0117 private:
0118 G4AllITFinder();
0119 static G4ThreadLocal G4AllITFinder* fpInstance;
0120 std::map<G4ITType, G4VITFinder*> fITSubManager;
0121
0122 int fVerbose;
0123 };
0124
0125 template<typename T>
0126 inline G4KDTreeResultHandle G4AllITFinder::FindNearest(const G4ThreeVector& pos,
0127 const T* it)
0128 {
0129 return G4ITFinder<T>::Instance()->FindNearest(pos, it);
0130 }
0131
0132 template<typename T>
0133 inline G4KDTreeResultHandle G4AllITFinder::FindNearest(const T* it0,
0134 const T* it)
0135 {
0136 return G4ITFinder<T>::Instance()->FindNearest(it0, it);
0137 }
0138
0139 template<typename T>
0140 inline G4KDTreeResultHandle G4AllITFinder::FindNearestInRange(const G4ThreeVector& pos,
0141 const T* it,
0142 G4double range)
0143 {
0144 return G4ITFinder<T>::Instance()->FindNearestInRange(pos, it, range);
0145 }
0146
0147 template<typename T>
0148 inline G4KDTreeResultHandle G4AllITFinder::FindNearestInRange(const T* it0,
0149 const T* it,
0150 G4double range)
0151 {
0152 return G4ITFinder<T>::Instance()->FindNearestInRange(it0, it, range);
0153 }
0154
0155 #endif