Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4KDTree.hh
Go to the documentation of this file.
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// $Id: G4KDTree.hh 64057 2012-10-30 15:04:49Z gcosmo $
27//
28// Author: Mathieu Karamitros (kara (AT) cenbg . in2p3 . fr)
29//
30// WARNING : This class is released as a prototype.
31// It might strongly evolve or even disapear in the next releases.
32//
33// History:
34// -----------
35// 10 Oct 2011 M.Karamitros created
36//
37// -------------------------------------------------------------------
38
39#ifndef G4KDTREE_HH
40#define G4KDTREE_HH
41
42#include <vector>
43#include "G4KDTreeResult.hh"
44
45//__________________________________
46// Methods to act on kdnode
47// Methods defined in G4KDNode.cc :
49void Free(G4KDNode*&);
50void* GetData(G4KDNode*);
51const double* GetNodePosition(G4KDNode*);
52//__________________________________
53
54/**
55 * G4KDTree is used by the ITManager to locate the neareast neighbours.
56 * A kdtree sorts out node in such a way that it reduces the number of node check.
57 * The results of this search can be retrieved by G4KDTreeResultHandle.
58 */
59
61{
62 friend class G4KDNode ;
63 int fDim;
64 struct HyperRect *fRect;
65 void (*fDestr)(void*);
66 int fNbNodes;
67
68protected :
70
71public :
72 G4KDTree(int dim = 3);
73 virtual ~G4KDTree();
74
75 void Clear();
76
77 inline int GetDim();
78 inline void SetDataDestructor(void (*fDestr)(void*));
79
80 int GetNbNodes() { return fNbNodes; }
81 G4KDNode* GetRoot() { return fRoot ; }
82
83 // Insert and attache the data to a node at the specified position
84 // In return, it gives you the corresponding node
85 G4KDNode* Insert(const double *pos, void *data);
86 G4KDNode* Insert(const double& x, const double& y, const double& z, void *data); // 3D
87
88 /* Find one of the nearest nodes from the specified point.
89 *
90 * This function returns a pointer to a result set with at most one element.
91 */
92 G4KDTreeResultHandle Nearest( const double *pos);
93 G4KDTreeResultHandle Nearest( const double& x, const double& y, const double& z); // 3D
95
96 /* Find any nearest nodes from the specified point within a range.
97 *
98 * This function returns a pointer to a result set, which can be manipulated
99 * by the G4KDTreeResult.
100 * The returned pointer can be null as an indication of an error. Otherwise
101 * a valid result set is always returned which may contain 0 or more elements.
102 */
103 G4KDTreeResultHandle NearestInRange( const double *pos, const double& range);
104 G4KDTreeResultHandle NearestInRange( const double& x,
105 const double& y,
106 const double& z,
107 const double& range); // 3D
108 G4KDTreeResultHandle NearestInRange( G4KDNode* node, const double& range);
109
110protected :
111 void __Clear_Rec(G4KDNode *node) ;
112
113 int __NearestInRange(G4KDNode *node,
114 const double *pos,
115 const double& range_sq,
116 const double& range,
117 G4KDTreeResult& list,
118 int ordered,
119 G4KDNode *source_node = 0);
120
121 void __NearestToPosition(G4KDNode *node,
122 const double *pos,
123 G4KDNode *&result,
124 double *result_dist_sq,
125 struct HyperRect* fRect);
126
127 void __NearestToNode(G4KDNode *source_node,
128 G4KDNode *node,
129 const double *pos,
130 std::vector<G4KDNode*>& result,
131 double *result_dist_sq,
132 struct HyperRect* fRect,
133 int& nbresult) ;
134};
135
137{
138 return fDim ;
139}
140
141void G4KDTree::SetDataDestructor(void (*fct)(void*))
142{
143 fDestr = fct;
144}
145
146#endif // G4KDTREE_HH
void * GetData(G4KDNode *)
Definition: G4KDNode.cc:45
const double * GetNodePosition(G4KDNode *)
Definition: G4KDNode.cc:50
void InactiveNode(G4KDNode *)
Definition: G4KDNode.cc:57
void Free(G4KDNode *&)
Definition: G4KDNode.cc:63
void __NearestToNode(G4KDNode *source_node, G4KDNode *node, const double *pos, std::vector< G4KDNode * > &result, double *result_dist_sq, struct HyperRect *fRect, int &nbresult)
Definition: G4KDTree.cc:420
G4KDNode * GetRoot()
Definition: G4KDTree.hh:81
void SetDataDestructor(void(*fDestr)(void *))
Definition: G4KDTree.hh:141
int GetDim()
Definition: G4KDTree.hh:136
int __NearestInRange(G4KDNode *node, const double *pos, const double &range_sq, const double &range, G4KDTreeResult &list, int ordered, G4KDNode *source_node=0)
Definition: G4KDTree.cc:261
G4KDNode * fRoot
Definition: G4KDTree.hh:69
void __NearestToPosition(G4KDNode *node, const double *pos, G4KDNode *&result, double *result_dist_sq, struct HyperRect *fRect)
Definition: G4KDTree.cc:308
G4KDTreeResultHandle NearestInRange(const double *pos, const double &range)
Definition: G4KDTree.cc:561
G4KDNode * Insert(const double *pos, void *data)
Definition: G4KDTree.cc:221
G4KDTreeResultHandle Nearest(const double *pos)
Definition: G4KDTree.cc:386
virtual ~G4KDTree()
Definition: G4KDTree.cc:178
int GetNbNodes()
Definition: G4KDTree.hh:80
void __Clear_Rec(G4KDNode *node)
Definition: G4KDTree.cc:203
void Clear()
Definition: G4KDTree.cc:190