Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4KDTree.cc
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.cc 64057 2012-10-30 15:04:49Z gcosmo $
27//
28// Author: Mathieu Karamitros (kara (AT) cenbg . in2p3 . fr)
29//
30// History:
31// -----------
32// 10 Oct 2011 M.Karamitros created
33//
34// -------------------------------------------------------------------
35/*
36 * Based on ``kdtree'', a library for working with kd-trees.
37 * Copyright (C) 2007-2009 John Tsiombikas <[email protected]>
38 * The original open-source version of this code
39 * may be found at http://code.google.com/p/kdtree/
40 *
41 * Redistribution and use in source and binary forms, with or without
42 * modification, are permitted provided that the following conditions are
43 * met:
44 * 1. Redistributions of source code must retain the above copyright
45 * notice, this
46 * list of conditions and the following disclaimer.
47 * 2. Redistributions in binary form must reproduce the above copyright
48 * notice,
49 * this list of conditions and the following disclaimer in the
50 * documentation
51 * and/or other materials provided with the distribution.
52 * 3. The name of the author may not be used to endorse or promote products
53 * derived from this software without specific prior written permission.
54 *
55 * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
56 * AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED
57 * WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE DISCLAIMED.
58 * IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT,
59 * INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT
60 * LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS;
61 * OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT,
62 * STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
63 * OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
64*/
65/* single nearest neighbor search written by Tamas Nepusz
67*/
68
69#include "globals.hh"
70#include <stdio.h>
71#include <stdlib.h>
72#include <cmath>
73#include "G4KDTree.hh"
74#include "G4KDNode.hh"
75#include "G4KDTreeResult.hh"
76#include <list>
77#include <iostream>
78
79using namespace std;
80
81//______________________________________________________________________
83{
84public:
85 HyperRect(int dim, const double *min, const double *max)
86 {
87 fDim = dim;
88 fMin = new double[fDim];
89 fMax = new double[fDim];
90 size_t size = fDim * sizeof(double);
91 memcpy(fMin, min, size);
92 memcpy(fMax, max, size);
93 }
94
95
97 {
98 delete[] fMin;
99 delete[] fMax;
100 }
101
102 HyperRect(const HyperRect& rect)
103 {
104 fDim = rect.fDim;
105 fMin = new double[fDim];
106 fMax = new double[fDim];
107 size_t size = fDim * sizeof(double);
108 memcpy(fMin, rect.fMin, size);
109 memcpy(fMax, rect.fMax, size);
110 }
111
112 void Extend(const double *pos)
113 {
114 int i;
115
116 for (i=0; i < fDim; i++)
117 {
118 if (pos[i] < fMin[i])
119 {
120 fMin[i] = pos[i];
121 }
122 if (pos[i] > fMax[i])
123 {
124 fMax[i] = pos[i];
125 }
126 }
127 }
128
129 bool CompareDistSqr(const double *pos, const double* bestmatch)
130 {
131 double result = 0;
132
133 for (int i=0; i < fDim; i++)
134 {
135 if (pos[i] < fMin[i])
136 {
137 result += sqr(fMin[i] - pos[i]);
138 }
139 else if (pos[i] > fMax[i])
140 {
141 result += sqr(fMax[i] - pos[i]);
142 }
143
144 if(result >= *bestmatch) return false ;
145 }
146
147 return true ;
148 }
149
150 int GetDim(){return fDim;}
151 double* GetMin(){return fMin;}
152 double* GetMax(){return fMax;}
153
154protected:
155 int fDim;
156 double *fMin, *fMax; /* minimum/maximum coords */
157
158private:
159 // should not be used
160 HyperRect& operator=(const HyperRect& rhs)
161 {
162 if(this == &rhs) return *this;
163 return *this;
164 }
165};
166
167//______________________________________________________________________
168// KDTree methods
170{
171 fDim = k;
172 fRoot = 0;
173 fDestr = 0;
174 fRect = 0;
175 fNbNodes = 0;
176}
177
179{
181 fRoot = 0;
182
183 if (fRect)
184 {
185 delete fRect;
186 fRect = 0;
187 }
188}
189
191{
193 fRoot = 0;
194 fNbNodes = 0;
195
196 if (fRect)
197 {
198 delete fRect;
199 fRect = 0;
200 }
201}
202
204{
205 if(!node) return;
206
207 if(node->GetLeft()) __Clear_Rec(node->GetLeft());
208 if(node->GetRight()) __Clear_Rec(node->GetRight());
209
210 if(fDestr)
211 {
212 if(node->GetData())
213 {
214 fDestr(node->GetData());
215 node->SetData(0);
216 }
217 }
218 delete node;
219}
220
221G4KDNode* G4KDTree::Insert(const double *pos, void *data)
222{
223 G4KDNode* node = 0 ;
224 if(!fRoot)
225 {
226 fRoot = new G4KDNode(this,pos,data,0, 0);
227 node = fRoot;
228 fNbNodes = 0;
229 fNbNodes++;
230 }
231 else
232 {
233 if((node=fRoot->Insert(pos, data)))
234 {
235 fNbNodes++;
236 }
237 }
238
239 if (fRect == 0)
240 {
241 fRect = new HyperRect(fDim,pos,pos);
242 }
243 else
244 {
245 fRect->Extend(pos);
246 }
247
248 return node;
249}
250
251G4KDNode* G4KDTree::Insert(const double& x, const double& y, const double& z, void *data)
252{
253 double buf[3];
254 buf[0] = x;
255 buf[1] = y;
256 buf[2] = z;
257 return Insert(buf, data);
258}
259
260//__________________________________________________________________
261int G4KDTree::__NearestInRange(G4KDNode *node, const double *pos, const double& range_sq,
262 const double& range, G4KDTreeResult& list, int ordered, G4KDNode *source_node)
263{
264 if(!node) return 0;
265
266 double dist_sq(DBL_MAX), dx(DBL_MAX);
267 int ret(-1), added_res(0);
268
269 if(node-> GetData() && node != source_node)
270 {
271 bool do_break = false ;
272 dist_sq = 0;
273 for(int i=0; i<fDim; i++)
274 {
275 dist_sq += sqr(node->GetPosition()[i] - pos[i]);
276 if(dist_sq > range_sq)
277 {
278 do_break = true;
279 break;
280 }
281 }
282 if(!do_break && dist_sq <= range_sq)
283 {
284 list.Insert(dist_sq, node);
285 added_res = 1;
286 }
287 }
288
289 dx = pos[node->GetAxis()] - node->GetPosition()[node->GetAxis()];
290
291 ret = __NearestInRange(dx <= 0.0 ? node->GetLeft() : node->GetRight(), pos, range_sq, range, list, ordered, source_node);
292 if(ret >= 0 && fabs(dx) <= range)
293 {
294 added_res += ret;
295 ret = __NearestInRange(dx <= 0.0 ? node->GetRight() : node->GetLeft(), pos, range_sq, range, list, ordered, source_node);
296 }
297
298 if(ret == -1)
299 {
300 return -1;
301 }
302 added_res += ret;
303
304 return added_res;
305}
306
307//__________________________________________________________________
308void G4KDTree::__NearestToPosition(G4KDNode *node, const double *pos, G4KDNode *&result,
309 double *result_dist_sq, HyperRect* rect)
310{
311 int dir = node->GetAxis();
312 int i;
313 double dummy(0.), dist_sq(-1.);
314 G4KDNode *nearer_subtree(0), *farther_subtree (0);
315 double *nearer_hyperrect_coord(0),*farther_hyperrect_coord(0);
316
317 /* Decide whether to go left or right in the tree */
318 dummy = pos[dir] - node->GetPosition()[dir];
319 if (dummy <= 0)
320 {
321 nearer_subtree = node->GetLeft();
322 farther_subtree = node->GetRight();
323
324 nearer_hyperrect_coord = rect->GetMax() + dir;
325 farther_hyperrect_coord = rect->GetMin() + dir;
326 }
327 else
328 {
329 nearer_subtree = node->GetRight();
330 farther_subtree = node->GetLeft();
331 nearer_hyperrect_coord = rect->GetMin() + dir;
332 farther_hyperrect_coord = rect->GetMax() + dir;
333 }
334
335 if (nearer_subtree)
336 {
337 /* Slice the hyperrect to get the hyperrect of the nearer subtree */
338 dummy = *nearer_hyperrect_coord;
339 *nearer_hyperrect_coord = node->GetPosition()[dir];
340 /* Recurse down into nearer subtree */
341 __NearestToPosition(nearer_subtree, pos, result, result_dist_sq, rect);
342 /* Undo the slice */
343 *nearer_hyperrect_coord = dummy;
344 }
345
346 /* Check the distance of the point at the current node, compare it
347 * with our best so far */
348 if(node->GetData())
349 {
350 dist_sq = 0;
351 bool do_break = false ;
352 for(i=0; i < fDim; i++)
353 {
354 dist_sq += sqr(node->GetPosition()[i] - pos[i]);
355 if(dist_sq > *result_dist_sq)
356 {
357 do_break = true;
358 break ;
359 }
360 }
361 if (!do_break && dist_sq < *result_dist_sq)
362 {
363 result = node;
364 *result_dist_sq = dist_sq;
365 }
366 }
367
368 if (farther_subtree)
369 {
370 /* Get the hyperrect of the farther subtree */
371 dummy = *farther_hyperrect_coord;
372 *farther_hyperrect_coord = node->GetPosition()[dir];
373 /* Check if we have to recurse down by calculating the closest
374 * point of the hyperrect and see if it's closer than our
375 * minimum distance in result_dist_sq. */
376 if (rect->CompareDistSqr(pos,result_dist_sq))
377 {
378 /* Recurse down into farther subtree */
379 __NearestToPosition(farther_subtree, pos, result, result_dist_sq, rect);
380 }
381 /* Undo the slice on the hyperrect */
382 *farther_hyperrect_coord = dummy;
383 }
384}
385
387{
388// G4cout << "Nearest(pos)" << G4endl ;
389
390 if (!fRect) return 0;
391
392 G4KDNode *result(0);
393 double dist_sq = DBL_MAX;
394
395 /* Duplicate the bounding hyperrectangle, we will work on the copy */
396 HyperRect *newrect = new HyperRect(*fRect);
397
398 /* Our first guesstimate is the root node */
399 /* Search for the nearest neighbour recursively */
400 __NearestToPosition(fRoot, pos, result, &dist_sq, newrect);
401
402 /* Free the copy of the hyperrect */
403 delete newrect;
404
405 /* Store the result */
406 if (result)
407 {
408 G4KDTreeResultHandle rset = new G4KDTreeResult(this);
409 rset->Insert(dist_sq, result);
410 rset -> Rewind();
411 return rset;
412 }
413 else
414 {
415 return 0;
416 }
417}
418
419//__________________________________________________________________
421 const double *pos, std::vector<G4KDNode*>& result, double *result_dist_sq,
422 HyperRect* rect, int& nbresult)
423{
424 int dir = node->GetAxis();
425 double dummy, dist_sq;
426 G4KDNode *nearer_subtree (0), *farther_subtree (0);
427 double *nearer_hyperrect_coord (0), *farther_hyperrect_coord (0);
428
429 /* Decide whether to go left or right in the tree */
430 dummy = pos[dir] - node->GetPosition()[dir];
431 if (dummy <= 0)
432 {
433 nearer_subtree = node->GetLeft();
434 farther_subtree = node->GetRight();
435 nearer_hyperrect_coord = rect->GetMax() + dir;
436 farther_hyperrect_coord = rect->GetMin() + dir;
437 }
438 else
439 {
440 nearer_subtree = node->GetRight();
441 farther_subtree = node->GetLeft();
442 nearer_hyperrect_coord = rect->GetMin() + dir;
443 farther_hyperrect_coord = rect->GetMax() + dir;
444 }
445
446 if (nearer_subtree)
447 {
448 /* Slice the hyperrect to get the hyperrect of the nearer subtree */
449 dummy = *nearer_hyperrect_coord;
450 *nearer_hyperrect_coord = node->GetPosition()[dir];
451 /* Recurse down into nearer subtree */
452 __NearestToNode(source_node, nearer_subtree, pos, result, result_dist_sq, rect, nbresult);
453 /* Undo the slice */
454 *nearer_hyperrect_coord = dummy;
455 }
456
457 /* Check the distance of the point at the current node, compare it
458 * with our best so far */
459 if(node->GetData() && node != source_node)
460 {
461 dist_sq = 0;
462 bool do_break = false;
463 for(int i=0; i < fDim; i++)
464 {
465 dist_sq += sqr(node->GetPosition()[i] - pos[i]);
466 if(dist_sq > *result_dist_sq)
467 {
468 do_break = true;
469 break ;
470 }
471 }
472 if(!do_break)
473 {
474 if (dist_sq < *result_dist_sq)
475 {
476 result.clear();
477 nbresult = 1 ;
478 result.push_back(node);
479 *result_dist_sq = dist_sq;
480 }
481 else if(dist_sq == *result_dist_sq)
482 {
483 result.push_back(node);
484 nbresult++;
485 }
486 }
487 }
488
489 if (farther_subtree)
490 {
491 /* Get the hyperrect of the farther subtree */
492 dummy = *farther_hyperrect_coord;
493 *farther_hyperrect_coord = node->GetPosition()[dir];
494 /* Check if we have to recurse down by calculating the closest
495 * point of the hyperrect and see if it's closer than our
496 * minimum distance in result_dist_sq. */
497 // if (hyperrect_dist_sq(rect, pos) < *result_dist_sq)
498 if (rect->CompareDistSqr(pos,result_dist_sq))
499 {
500 /* Recurse down into farther subtree */
501 __NearestToNode(source_node, farther_subtree, pos, result, result_dist_sq, rect, nbresult);
502 }
503 /* Undo the slice on the hyperrect */
504 *farther_hyperrect_coord = dummy;
505 }
506}
507
509{
510// G4cout << "Nearest(node)" << G4endl ;
511 if (!fRect)
512 {
513 G4cout << "Tree empty" << G4endl ;
514 return 0;
515 }
516
517 const double* pos = node->GetPosition();
518 std::vector<G4KDNode*> result;
519 double dist_sq = DBL_MAX;
520
521 /* Duplicate the bounding hyperrectangle, we will work on the copy */
522 HyperRect *newrect = new HyperRect(*fRect);
523
524 /* Search for the nearest neighbour recursively */
525 int nbresult = 0 ;
526
527 __NearestToNode(node, fRoot, pos, result, &dist_sq, newrect, nbresult);
528
529 /* Free the copy of the hyperrect */
530 delete newrect;
531
532 /* Store the result */
533 if (!result.empty())
534 {
535 G4KDTreeResultHandle rset(new G4KDTreeResult(this));
536 int j = 0 ;
537 while (j<nbresult)
538 {
539 rset->Insert(dist_sq, result[j]);
540 j++;
541 }
542 rset->Rewind();
543
544 return rset;
545 }
546 else
547 {
548 return 0;
549 }
550}
551
552G4KDTreeResultHandle G4KDTree::Nearest(const double& x, const double& y, const double& z)
553{
554 double pos[3];
555 pos[0] = x;
556 pos[1] = y;
557 pos[2] = z;
558 return Nearest(pos);
559}
560
561G4KDTreeResultHandle G4KDTree::NearestInRange(const double *pos, const double& range)
562{
563 int ret(-1);
564
565 const double range_sq = sqr(range) ;
566
567 G4KDTreeResultHandle rset = new G4KDTreeResult(this);
568 if((ret = __NearestInRange(fRoot, pos, range_sq, range, *(rset()), 0)) == -1)
569 {
570 rset = 0;
571 return rset;
572 }
573 rset->Sort();
574 rset->Rewind();
575 return rset;
576}
577
579 const double& y,
580 const double& z,
581 const double& range)
582{
583 double buf[3];
584 buf[0] = x;
585 buf[1] = y;
586 buf[2] = z;
587 return NearestInRange(buf, range);
588}
589
591{
592 if(!node) return 0 ;
593 int ret(-1);
594
595 G4KDTreeResult *rset = new G4KDTreeResult(this);
596
597 const double range_sq = sqr(range) ;
598
599 if((ret = __NearestInRange(fRoot, node->GetPosition(), range_sq, range, *rset, 0, node)) == -1)
600 {
601 delete rset;
602 return 0;
603 }
604 rset->Sort();
605 rset->Rewind();
606 return rset;
607}
void * GetData(G4KDNode *)
Definition: G4KDNode.cc:45
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
G4KDNode * GetLeft()
Definition: G4KDNode.hh:137
void SetData(void *)
Definition: G4KDNode.hh:122
G4KDNode * GetRight()
Definition: G4KDNode.hh:142
const double * GetPosition()
Definition: G4KDNode.hh:127
int GetAxis()
Definition: G4KDNode.hh:112
void * GetData()
Definition: G4KDNode.hh:117
G4KDNode * Insert(const double *p, void *data)
Definition: G4KDNode.cc:156
void Insert(double, G4KDNode *)
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
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
friend class G4KDNode
Definition: G4KDTree.hh:62
G4KDTreeResultHandle NearestInRange(const double *pos, const double &range)
Definition: G4KDTree.cc:561
G4KDTree(int dim=3)
Definition: G4KDTree.cc:169
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
void __Clear_Rec(G4KDNode *node)
Definition: G4KDTree.cc:203
void Clear()
Definition: G4KDTree.cc:190
int GetDim()
Definition: G4KDTree.cc:150
~HyperRect()
Definition: G4KDTree.cc:96
HyperRect(const HyperRect &rect)
Definition: G4KDTree.cc:102
int fDim
Definition: G4KDTree.cc:155
double * fMin
Definition: G4KDTree.cc:156
bool CompareDistSqr(const double *pos, const double *bestmatch)
Definition: G4KDTree.cc:129
double * GetMax()
Definition: G4KDTree.cc:152
double * fMax
Definition: G4KDTree.cc:156
void Extend(const double *pos)
Definition: G4KDTree.cc:112
double * GetMin()
Definition: G4KDTree.cc:151
HyperRect(int dim, const double *min, const double *max)
Definition: G4KDTree.cc:85
T sqr(const T &x)
Definition: templates.hh:145
#define DBL_MAX
Definition: templates.hh:83