FastJet 3.5.0
Loading...
Searching...
No Matches
DnnPlane.hh
1//FJSTARTHEADER
2// $Id$
3//
4// Copyright (c) 2005-2025, Matteo Cacciari, Gavin P. Salam and Gregory Soyez
5//
6//----------------------------------------------------------------------
7// This file is part of FastJet.
8//
9// FastJet is free software; you can redistribute it and/or modify
10// it under the terms of the GNU General Public License as published by
11// the Free Software Foundation; either version 2 of the License, or
12// (at your option) any later version.
13//
14// The algorithms that underlie FastJet have required considerable
15// development. They are described in the original FastJet paper,
16// hep-ph/0512210 and in the manual, arXiv:1111.6097. If you use
17// FastJet as part of work towards a scientific publication, please
18// quote the version you use and include a citation to the manual and
19// optionally also to hep-ph/0512210.
20//
21// FastJet is distributed in the hope that it will be useful,
22// but WITHOUT ANY WARRANTY; without even the implied warranty of
23// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
24// GNU General Public License for more details.
25//
26// You should have received a copy of the GNU General Public License
27// along with FastJet. If not, see <http://www.gnu.org/licenses/>.
28//----------------------------------------------------------------------
29//FJENDHEADER
30
31
32#include "fastjet/config.h"
33
34#ifndef DROP_CGAL // in case we do not have the code for CGAL
35
36#ifndef __FASTJET_DNNPLANE_HH__
37#define __FASTJET_DNNPLANE_HH__
38
39#include "fastjet/internal/Triangulation.hh"
40#include "fastjet/internal/DynamicNearestNeighbours.hh"
41
42FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh
43
44
45/// \if internal_doc
46/// @ingroup internal
47/// \class DnnPlane
48/// class derived from DynamicNearestNeighbours that provides an
49/// implementation for the Euclidean plane
50///
51/// This class that uses CGAL Delaunay triangulation for most of the
52/// work (it allows for easy and efficient removal and addition of
53/// points and circulation over a point's neighbours). The treatment
54/// of coincident points is not supported by CGAL and is implemented
55/// according to the method specified in
56/// issue-tracker/2012-02-CGAL-coincident/METHOD
57/// \endif
59 public:
60 /// empty initaliser
62
63 /// Initialiser from a set of points on an Eta-Phi plane, where both
64 /// eta and phi can have arbitrary ranges
65 DnnPlane(const std::vector<EtaPhi> &, const bool & verbose = false );
66
67
68 /// Returns the index of the nearest neighbour of point labelled
69 /// by ii (assumes ii is valid)
70 int NearestNeighbourIndex(const int ii) const ;
71
72 /// Returns the distance to the nearest neighbour of point labelled
73 /// by index ii (assumes ii is valid)
74 double NearestNeighbourDistance(const int ii) const ;
75
76 /// Returns true iff the given index corresponds to a point that
77 /// exists in the DNN structure (meaning that it has been added, and
78 /// not removed in the meantime)
79 bool Valid(const int index) const;
80
81 void RemoveAndAddPoints(const std::vector<int> & indices_to_remove,
82 const std::vector<EtaPhi> & points_to_add,
83 std::vector<int> & indices_added,
84 std::vector<int> & indices_of_updated_neighbours);
85
86 /// returns the EtaPhi of point with index i.
87 EtaPhi etaphi(const int i) const;
88 /// returns the eta point with index i.
89 double eta(const int i) const;
90 /// returns the phi point with index i.
91 double phi(const int i) const;
92
93private:
94
95 /// Structure containing a vertex_handle and cached information on
96 /// the nearest neighbour.
97 struct SuperVertex {
98 Vertex_handle vertex; // NULL indicates inexistence...
99 double NNdistance;
100 int NNindex;
101 int coincidence; // ==vertex->info.val() if no coincidence
102 // points to the coinciding SV in case of coincidence
103 // later on for cylinder put a second vertex?
104 };
105
106 std::vector<SuperVertex> _supervertex;
107 //set<Vertex_handle> _vertex_set;
108 bool _verbose;
109
110 //static const bool _crash_on_coincidence = true;
111 static const bool _crash_on_coincidence = false;
112
113 Triangulation _TR; /// CGAL object for dealing with triangulations
114
115 /// calculates and returns the euclidean distance between points p1
116 /// and p2
117 inline double _euclid_distance(const Point& p1, const Point& p2) const {
118 double distx= p1.x()-p2.x();
119 double disty= p1.y()-p2.y();
120 return distx*distx+disty*disty;
121 }
122
123 //----------------------------------------------------------------------
124 /// Determines the index and distance of the nearest neighbour to
125 /// point j and puts the information into the _supervertex entry for j
126 void _SetNearest(const int j);
127
128 //----------------------------------------------------------------------
129 /// Determines and stores the nearest neighbour of j.
130 ///
131 /// For each voronoi neighbour D of j if the distance between j and D
132 /// is less than D's own nearest neighbour, then update the
133 /// nearest-neighbour info in D; push D's index onto
134 /// indices_of_updated_neighbours
135 ///
136 /// Note that j is NOT pushed onto indices_of_updated_neighbours --
137 /// if you want it there, put it there yourself.
138 void _SetAndUpdateNearest(const int j,
139 std::vector<int> & indices_of_updated_neighbours);
140
141 /// given a vertex_handle returned by CGAL on insertion of a new
142 /// points, returns the coinciding vertex's value if it turns out
143 /// that it corresponds to a vertex that we already knew about
144 /// (usually because two points coincide)
145 int _CheckIfVertexPresent(const Vertex_handle & vertex,
146 const int its_index);
147
148 //----------------------------------------------------------------------
149 /// if the distance between 'pref' and 'candidate' is smaller (or
150 /// equal) than the one between 'pref' and 'near', return true and
151 /// set 'mindist' to that distance. Note that it is assumed that
152 /// 'mindist' is the euclidian distance between 'pref' and 'near'
153 ///
154 /// Note that the 'near' point is passed through its vertex rather
155 /// than as a point. This allows us to handle cases where we have no min
156 /// yet (near is the infinite vertex)
157 inline bool _is_closer_to(const Point &pref,
158 const Point &candidate,
159 const Vertex_handle &near,
160 double & dist,
161 double & mindist){
162 dist = _euclid_distance(pref, candidate);
163 return _is_closer_to_with_hint(pref, candidate, near, dist, mindist);
164 }
165
166 /// same as '_is_closer_to' except that 'dist' already contains the
167 /// distance between 'pref' and 'candidate'
168 inline bool _is_closer_to_with_hint(const Point &pref,
169 const Point &candidate,
170 const Vertex_handle &near,
171 const double & dist,
172 double & mindist){
173
174 // check if 'dist', the pre-computed distance between 'candidate'
175 // and 'pref' is smaller than the distance between 'pref' and its
176 // currently registered nearest neighbour 'near' (and update
177 // things if it is)
178 //
179 // Interestingly enough, it has to be pointed out that the use of
180 // 'abs' instead of 'std::abs' returns wrong results (apparently
181 // ints without any compiler warning)
182 //
183 // The (near != NULL) test is there for one single reason: when
184 // checking that a newly inserted point is not closer than a
185 // previous NN, if that distance comparison involves a "nearly
186 // degenerate" distance we need to access near->point. But
187 // sometimes, in the course of RemoveAndAddPoints, its previous NN
188 // has been deleted and its vertex (corresponding to 'near') set
189 // to NULL. This is not a problem as all points having a deleted
190 // point as NN will have their NN explicitly recomputed at the end
191 // of RemoveAndAddPoints so here we should just make sure there is
192 // no crash... that's done by checking (near != NULL)
193 if ((std::abs(dist-mindist)<DISTANCE_FOR_CGAL_CHECKS) &&
194 _is_not_null(near) && // GPS cf. note below about != NULL with clang/CGAL
195 (_euclid_distance(candidate, near->point())<DISTANCE_FOR_CGAL_CHECKS)){
196 // we're in a situation where there might be a rounding issue,
197 // use CGAL's distance computation to get it right
198 //
199 // Note that in the test right above,
200 // (abs(dist-mindist)<1e-12) guarantees that the current
201 // nearest point is not the infinite vertex and thus
202 // nearest->point() is not ill-defined
203 if (_verbose) std::cout << "using CGAL's distance ordering" << std::endl;
204 if (CGAL::compare_distance_to_point(pref, candidate, near->point())!=CGAL::LARGER){
205 mindist = dist;
206 return true;
207 }
208 } else if (dist <= mindist) {
209 // Note that the use of a <= in the above expression (instead of
210 // a strict ordering <) is important in one case: when checking
211 // if a new point is the new NN of one of the points in its
212 // neighbourhood, in case of distances being ==, we are sure
213 // that 'candidate' is in a cell adjacent to 'pref' while it may
214 // no longer be the case for 'near'
215 mindist = dist;
216 return true;
217 }
218
219 return false;
220 }
221
222 /// if a distance between a point and 2 others is smaller than this
223 /// and the distance between the two points is also smaller than this
224 /// then use CGAL to compare the distances.
225 static const double DISTANCE_FOR_CGAL_CHECKS;
226
227
228 /// small routine to check if if a vertex handle is not null.
229 ///
230 /// This is centralised here because of the need for a not-very
231 /// readable workaround for certain CGAL/clang++ compiler
232 /// combinations.
233 inline static bool _is_not_null(const Vertex_handle & vh) {
234 // as found in issue 2015-07-clang61+CGAL, the form
235 // vh != NULL
236 // appears to be broken with CGAL-3.6 and clang-3.6.0/.1
237 //
238 // The not very "readable"
239 // vh.operator->()
240 // does the job instead
241 return vh.operator->();
242 }
243
244};
245
246
247// here follow some inline implementations of the simpler of the
248// functions defined above
249
250inline int DnnPlane::NearestNeighbourIndex(const int ii) const {
251 return _supervertex[ii].NNindex;}
252
253inline double DnnPlane::NearestNeighbourDistance(const int ii) const {
254 return _supervertex[ii].NNdistance;}
255
256inline bool DnnPlane::Valid(const int index) const {
257 if (index >= 0 && index < static_cast<int>(_supervertex.size())) {
258 return _is_not_null(_supervertex[index].vertex);
259 }
260 else {
261 return false;
262 }
263}
264
265
266inline EtaPhi DnnPlane::etaphi(const int i) const {
267 Point * p = & (_supervertex[i].vertex->point());
268 return EtaPhi(p->x(),p->y()); }
269
270inline double DnnPlane::eta(const int i) const {
271 return _supervertex[i].vertex->point().x(); }
272
273inline double DnnPlane::phi(const int i) const {
274 return _supervertex[i].vertex->point().y(); }
275
276
277FASTJET_END_NAMESPACE
278
279#endif // __FASTJET_DNNPLANE_HH__
280
281#endif // DROP_CGAL
DnnPlane(const std::vector< EtaPhi > &, const bool &verbose=false)
Initialiser from a set of points on an Eta-Phi plane, where both eta and phi can have arbitrary range...
DnnPlane()
empty initaliser
Definition DnnPlane.hh:61
Shortcut for dealing with eta-phi coordinates.