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 // Author: Mathieu Karamitros 28 29 // The code is developed in the framework of the ESA AO7146 30 // 31 // We would be very happy hearing from you, send us your feedback! :) 32 // 33 // In order for Geant4-DNA to be maintained and still open-source, 34 // article citations are crucial. 35 // If you use Geant4-DNA chemistry and you publish papers about your software, 36 // in addition to the general paper on Geant4-DNA: 37 // 38 // Int. J. Model. Simul. Sci. Comput. 1 (2010) 157–178 39 // 40 // we would be very happy if you could please also cite the following 41 // reference papers on chemistry: 42 // 43 // J. Comput. Phys. 274 (2014) 841-882 44 // Prog. Nucl. Sci. Tec. 2 (2011) 503-508 45 46 #ifndef G4KDTREE_HH 47 #define G4KDTREE_HH 1 48 49 #include <vector> 50 #include "G4Types.hh" 51 #include "G4KDNode.hh" 52 #include "G4KDTreeResult.hh" 53 54 class G4KDMap; 55 template <typename PointT> 56 class G4KDNode; 57 58 //__________________________________ 59 // Methods to act on kdnode 60 // Methods defined in G4KDNode.cc : 61 void InactiveNode(G4KDNode_Base*); 62 void Free(G4KDNode_Base*&); 63 //__________________________________ 64 65 /** 66 * G4KDTree is used by the ITManager to locate the neareast neighbours. 67 * A kdtree sorts out node in such a way that it reduces the number of node 68 * check. The results of this search can be retrieved by G4KDTreeResultHandle. 69 */ 70 class G4KDTree 71 { 72 friend class G4KDNode_Base; 73 74 public: 75 G4KDTree(std::size_t dim = 3); 76 ~G4KDTree(); 77 void Clear(); 78 79 void Print(std::ostream& out = G4cout) const; 80 void Build(); 81 void NoticeNodeDeactivation() 82 { 83 fNbActiveNodes--; 84 if(fNbActiveNodes <= 0) 85 { 86 Clear(); 87 } 88 } 89 90 std::size_t GetDim() const { return fDim; } 91 G4int GetNbNodes() const { return fNbNodes; } 92 G4KDNode_Base* GetRoot() { return fRoot; } 93 94 template <typename PointT> 95 G4KDNode_Base* InsertMap(PointT* pos); 96 97 // Insert and attache the data to a node at the specified position 98 // In return, it gives you the corresponding node 99 template <typename PointT> 100 G4KDNode_Base* Insert(PointT* pos); // 3D 101 102 template <typename PointT> 103 G4KDNode_Base* Insert(const PointT& pos); // 3D 104 105 /* Find one of the nearest nodes from the specified point. 106 * 107 * This function returns a pointer to a result set with at most one element. 108 */ 109 template <typename Position> 110 G4KDTreeResultHandle Nearest(const Position& pos); 111 G4KDTreeResultHandle Nearest(G4KDNode_Base* node); 112 113 /* Find any nearest nodes from the specified point within a range. 114 * 115 * This function returns a pointer to a result set, which can be manipulated 116 * by the G4KDTreeResult. 117 * The returned pointer can be null as an indication of an error. Otherwise 118 * a valid result set is always returned which may contain 0 or more elements. 119 */ 120 template <typename Position> 121 G4KDTreeResultHandle NearestInRange(const Position& pos, const G4double& range); 122 G4KDTreeResultHandle NearestInRange(G4KDNode_Base* node, const G4double& range); 123 124 void* operator new(std::size_t); 125 void operator delete(void*); 126 127 protected: 128 //______________________________________________________________________ 129 class HyperRect 130 { 131 public: 132 HyperRect(std::size_t dim) 133 : fDim(dim) 134 , fMin(new G4double[fDim]) 135 , fMax(new G4double[fDim]) 136 {} 137 138 template <typename Position> 139 void SetMinMax(const Position& min, const Position& max) 140 { 141 for(std::size_t i = 0; i < fDim; ++i) 142 { 143 fMin[i] = min[(G4int)i]; 144 fMax[i] = max[(G4int)i]; 145 } 146 } 147 148 ~HyperRect() 149 { 150 delete[] fMin; 151 delete[] fMax; 152 } 153 154 HyperRect(const HyperRect& rect) 155 { 156 fDim = rect.fDim; 157 fMin = new G4double[fDim]; 158 fMax = new G4double[fDim]; 159 160 for(std::size_t i = 0; i < fDim; ++i) 161 { 162 fMin[i] = rect.fMin[i]; 163 fMax[i] = rect.fMax[i]; 164 } 165 } 166 167 template <typename Position> 168 void Extend(const Position& pos) 169 { 170 for(G4int i = 0; i < (G4int)fDim; ++i) 171 { 172 if(pos[i] < fMin[i]) 173 { 174 fMin[i] = pos[i]; 175 } 176 if(pos[i] > fMax[i]) 177 { 178 fMax[i] = pos[i]; 179 } 180 } 181 } 182 183 template <typename Position> 184 G4bool CompareDistSqr(const Position& pos, const G4double* bestmatch) 185 { 186 G4double result = 0; 187 188 for(std::size_t i = 0; i < fDim; ++i) 189 { 190 if(pos[(G4int)i] < fMin[i]) 191 { 192 result += sqr(fMin[i] - pos[(G4int)i]); 193 } 194 else if(pos[(G4int)i] > fMax[i]) 195 { 196 result += sqr(fMax[i] - pos[(G4int)i]); 197 } 198 199 if(result >= *bestmatch){ 200 return false; 201 } 202 } 203 204 return true; 205 } 206 207 std::size_t GetDim() { return fDim; } 208 G4double* GetMin() { return fMin; } 209 G4double* GetMax() { return fMax; } 210 211 protected: 212 std::size_t fDim; 213 G4double *fMin, *fMax; /* minimum/maximum coords */ 214 215 private: 216 // should not be used 217 HyperRect& operator=(const HyperRect& rhs) 218 { 219 if (this == &rhs) return *this; 220 return *this; 221 } 222 }; 223 224 protected: 225 void __InsertMap(G4KDNode_Base* node); 226 void __Clear_Rec(G4KDNode_Base* node); 227 228 template <typename Position> 229 G4int __NearestInRange(G4KDNode_Base* node, const Position& pos, 230 const G4double& range_sq, const G4double& range, 231 G4KDTreeResult& list, G4int ordered, 232 G4KDNode_Base* source_node = nullptr); 233 234 template <typename Position> 235 void __NearestToPosition(G4KDNode_Base* node, const Position& pos, 236 G4KDNode_Base*& result, G4double* result_dist_sq, 237 HyperRect* fRect); 238 239 template <typename Position> 240 void __NearestToNode(G4KDNode_Base* source_node, G4KDNode_Base* node, 241 const Position& pos, std::vector<G4KDNode_Base*>& result, 242 G4double* result_dist_sq, HyperRect* fRect, G4int& nbresult); 243 244 protected: 245 HyperRect* fRect = nullptr; 246 G4KDNode_Base* fRoot = nullptr; 247 std::size_t fDim; 248 G4int fNbNodes = 0; 249 G4int fNbActiveNodes = 0; 250 G4KDMap* fKDMap; 251 static G4Allocator<G4KDTree>*& fgAllocator(); 252 }; 253 254 #include "G4KDTree.icc" 255 256 #endif // G4KDTREE_HH 257