[11] | 1 | //STARTHEADER
|
---|
| 2 | // $Id: DnnPlane.cc,v 1.1 2008-11-06 14:32:14 ovyn Exp $
|
---|
| 3 | //
|
---|
| 4 | // Copyright (c) 2005-2006, Matteo Cacciari and Gavin Salam
|
---|
| 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 and are described in hep-ph/0512210. If you use
|
---|
| 16 | // FastJet as part of work towards a scientific publication, please
|
---|
| 17 | // include a citation to the FastJet paper.
|
---|
| 18 | //
|
---|
| 19 | // FastJet is distributed in the hope that it will be useful,
|
---|
| 20 | // but WITHOUT ANY WARRANTY; without even the implied warranty of
|
---|
| 21 | // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
---|
| 22 | // GNU General Public License for more details.
|
---|
| 23 | //
|
---|
| 24 | // You should have received a copy of the GNU General Public License
|
---|
| 25 | // along with FastJet; if not, write to the Free Software
|
---|
| 26 | // Foundation, Inc.:
|
---|
| 27 | // 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
|
---|
| 28 | //----------------------------------------------------------------------
|
---|
| 29 | //ENDHEADER
|
---|
| 30 |
|
---|
| 31 |
|
---|
| 32 | #ifndef DROP_CGAL // in case we do not have the code for CGAL
|
---|
| 33 |
|
---|
| 34 | #include<set>
|
---|
| 35 | #include<list>
|
---|
| 36 | #include "fastjet/internal/DnnPlane.hh"
|
---|
| 37 | using namespace std;
|
---|
| 38 |
|
---|
| 39 | FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh
|
---|
| 40 |
|
---|
| 41 |
|
---|
| 42 | /// Initialiser from a set of points on an Eta-Phi plane, where both
|
---|
| 43 | /// eta and phi can have arbitrary ranges
|
---|
| 44 | DnnPlane::DnnPlane(const vector<EtaPhi> & input_points,
|
---|
| 45 | const bool & verbose ) {
|
---|
| 46 |
|
---|
| 47 | _verbose = verbose;
|
---|
| 48 | int n = input_points.size();
|
---|
| 49 |
|
---|
| 50 | // construct Voronoi diagram in such a way as to get the vertex handles
|
---|
| 51 | // and remember to set CGAL info with the index of the vertex
|
---|
| 52 | SuperVertex sv;
|
---|
| 53 | for (int i = 0; i < n; i++) {
|
---|
| 54 | sv.vertex =
|
---|
| 55 | _TR.insert(Point(input_points[i].first, input_points[i].second));
|
---|
| 56 |
|
---|
| 57 | // we are not up to dealing with coincident vertices, so make
|
---|
| 58 | // sure the user knows!
|
---|
| 59 | _CrashIfVertexPresent(sv.vertex, i);
|
---|
| 60 |
|
---|
| 61 | // we need to assicate an index to each vertex -- thus when we get
|
---|
| 62 | // a vertex (e.g. as a nearest neighbour) from CGAL, we will be
|
---|
| 63 | // able to figure out which particle it corresponded to.
|
---|
| 64 | sv.vertex->info() = i;
|
---|
| 65 | _supervertex.push_back(sv);
|
---|
| 66 | }
|
---|
| 67 |
|
---|
| 68 | // label infinite vertex info with negative index
|
---|
| 69 | _TR.infinite_vertex()->info() = INFINITE_VERTEX;
|
---|
| 70 |
|
---|
| 71 | // set up the structure that holds nearest distances and neighbours
|
---|
| 72 | for (int j = 0; j < n; j++) {_SetNearest(j);}
|
---|
| 73 |
|
---|
| 74 | }
|
---|
| 75 |
|
---|
| 76 |
|
---|
| 77 | //----------------------------------------------------------------------
|
---|
| 78 | /// Crashes if the given vertex handle already exists. Otherwise
|
---|
| 79 | /// it does the bookkeeping for future such tests
|
---|
| 80 | void DnnPlane::_CrashIfVertexPresent(
|
---|
| 81 | const Vertex_handle & vertex, const int & its_index) {
|
---|
| 82 | if (!_crash_on_coincidence) return;
|
---|
| 83 |
|
---|
| 84 | // vertices that do not have the same geometric position as any
|
---|
| 85 | // other vertex so far added have info().val() == NEW_VERTEX -- this
|
---|
| 86 | // is ensured by the InitialisedInt class, which forms the "info"
|
---|
| 87 | // part of our
|
---|
| 88 | // CGAL::Triangulation_vertex_base_with_info_2<InitialisedInt,K>
|
---|
| 89 | //
|
---|
| 90 | // If the vertex coincides with one that already exists, then
|
---|
| 91 | // info().val() it's info().val() will have been updated (in
|
---|
| 92 | // DNN:DNN) to be equal to a vertex "index".
|
---|
| 93 | if (vertex->info().val() != NEW_VERTEX) {
|
---|
| 94 | ostringstream err;
|
---|
| 95 | err << "ERROR in DnnPlane::_CrashIfVertexPresent"
|
---|
| 96 | <<endl << "Point "<<its_index<<" coincides with point "
|
---|
| 97 | <<vertex->info().val() << endl;
|
---|
| 98 | throw DnnError(err.str());
|
---|
| 99 | }
|
---|
| 100 | }
|
---|
| 101 |
|
---|
| 102 |
|
---|
| 103 | //----------------------------------------------------------------------
|
---|
| 104 | /// remove the points labelled by the vector indices_to_remove, and
|
---|
| 105 | /// add the points specified by the vector points_to_add
|
---|
| 106 | /// (corresponding indices will be calculated automatically); the
|
---|
| 107 | /// idea behind this routine is that the points to be added will
|
---|
| 108 | /// somehow be close to the one or other of the points being removed
|
---|
| 109 | /// and this can be used by the implementation to provide hints for
|
---|
| 110 | /// inserting the new points in whatever structure it is using. In a
|
---|
| 111 | /// kt-algorithm the points being added will be a result of a
|
---|
| 112 | /// combination of the points to be removed -- hence the proximity
|
---|
| 113 | /// is (more or less) guaranteed.
|
---|
| 114 | void DnnPlane::RemoveAndAddPoints(
|
---|
| 115 | const vector<int> & indices_to_remove,
|
---|
| 116 | const vector<EtaPhi> & points_to_add,
|
---|
| 117 | vector<int> & indices_added,
|
---|
| 118 | vector<int> & indices_of_updated_neighbours) {
|
---|
| 119 |
|
---|
| 120 |
|
---|
| 121 | // build set of UNION of Voronoi neighbours of a pair of nearest
|
---|
| 122 | // neighbours
|
---|
| 123 | set<int> NeighbourUnion;
|
---|
| 124 | // later on it will be convenient to have access to a set (rather
|
---|
| 125 | // than vector) of indices being removed
|
---|
| 126 | set<int> indices_removed;
|
---|
| 127 |
|
---|
| 128 | // for each of the indices to be removed add the voronoi neighbourhood to
|
---|
| 129 | // the NeighbourUnion set.
|
---|
| 130 | for (size_t ir = 0; ir < indices_to_remove.size(); ir++) {
|
---|
| 131 | int index = indices_to_remove[ir];
|
---|
| 132 | indices_removed.insert(index);
|
---|
| 133 | if (_verbose) cout << " Starting RemoveAndAddPoints" << endl;
|
---|
| 134 | if (_verbose) cout << " point " << index << endl;
|
---|
| 135 | // have a circulators that will go round the Voronoi neighbours of
|
---|
| 136 | // _supervertex[index1].vertex
|
---|
| 137 | Vertex_circulator vc = _TR.incident_vertices(_supervertex[index].vertex);
|
---|
| 138 | Vertex_circulator done = vc;
|
---|
| 139 | do {
|
---|
| 140 | // if a neighbouring vertex not the infinite vertex, then add it
|
---|
| 141 | // to our union of neighbouring vertices.
|
---|
| 142 | if (_verbose) cout << "examining " << vc->info().val() << endl;
|
---|
| 143 | if (vc->info().val() != INFINITE_VERTEX) {
|
---|
| 144 | // NB: from it=1 onwards occasionally it might already have
|
---|
| 145 | // been inserted -- but double insertion still leaves only one
|
---|
| 146 | // copy in the set, so there's no problem
|
---|
| 147 | NeighbourUnion.insert(vc->info().val());
|
---|
| 148 | if (_verbose) cout << "inserted " << vc->info().val() << endl;
|
---|
| 149 | }
|
---|
| 150 | } while (++vc != done);
|
---|
| 151 | }
|
---|
| 152 |
|
---|
| 153 | if (_verbose) {
|
---|
| 154 | set<int>::iterator it = NeighbourUnion.begin();
|
---|
| 155 | cout << "Union of neighbours of combined points" << endl;
|
---|
| 156 | for ( ; it != NeighbourUnion.end(); ++it ) {
|
---|
| 157 | cout << *it << endl;
|
---|
| 158 | }
|
---|
| 159 | }
|
---|
| 160 |
|
---|
| 161 | // update set, triangulation and supervertex info
|
---|
| 162 | for (size_t ir = 0; ir < indices_to_remove.size(); ir++) {
|
---|
| 163 | int index = indices_to_remove[ir];
|
---|
| 164 |
|
---|
| 165 | // NeighbourUnion should not contain the points to be removed
|
---|
| 166 | // (because later we will assume they still exist).
|
---|
| 167 | NeighbourUnion.erase(indices_to_remove[ir]);
|
---|
| 168 |
|
---|
| 169 | // points to be removed should also be eliminated from the
|
---|
| 170 | // triangulation and the supervertex structure should be updated
|
---|
| 171 | // to reflect the fact that the points are no longer valid.
|
---|
| 172 | _TR.remove(_supervertex[index].vertex);
|
---|
| 173 | _supervertex[index].vertex = NULL;
|
---|
| 174 | }
|
---|
| 175 |
|
---|
| 176 | // add new point: give a "hint" to the inserter that
|
---|
| 177 | // the new point should be added close to old points -- the easiest way
|
---|
| 178 | // of getting this is to take a point from the NeighbourUnion AFTER we have
|
---|
| 179 | // removed point1, point2, and to get one of its incident faces.
|
---|
| 180 | //
|
---|
| 181 | // This hinting improves speed by c. 25% for 10^4 points because it
|
---|
| 182 | // avoids the need for a costly (sqrt{N}) location step (at least
|
---|
| 183 | // with a non-hierarchical triangulation -- with a hierarchical one,
|
---|
| 184 | // this step could be done away with, though there will still be a
|
---|
| 185 | // cost of O(ln N) to pay.
|
---|
| 186 | //
|
---|
| 187 | // For some reason inserting the point before the two removals
|
---|
| 188 | // slows things down by c. 25%. This importance of the order
|
---|
| 189 | // is not understood.
|
---|
| 190 | //
|
---|
| 191 | // At some point it might be worth trying to select the "nearest"
|
---|
| 192 | // of the various points in the neighbour union to avoid large
|
---|
| 193 | // steps in cases where we have 0..2pi periodicity and the first member
|
---|
| 194 | // of the neighbour union happens to be on the wrong side.
|
---|
| 195 | Face_handle face;
|
---|
| 196 | if (indices_to_remove.size() > 0) {
|
---|
| 197 | // face can only be found if there were points to remove in first place
|
---|
| 198 | face = _TR.incident_faces(
|
---|
| 199 | _supervertex[*NeighbourUnion.begin()].vertex);}
|
---|
| 200 | // make sure the output arrays are empty
|
---|
| 201 | indices_added.clear();
|
---|
| 202 | indices_of_updated_neighbours.clear();
|
---|
| 203 | for (size_t ia = 0; ia < points_to_add.size(); ia++) {
|
---|
| 204 | SuperVertex sv;
|
---|
| 205 | _supervertex.push_back(sv);
|
---|
| 206 | int index = _supervertex.size()-1;
|
---|
| 207 | indices_added.push_back(index);
|
---|
| 208 |
|
---|
| 209 | if (indices_to_remove.size() > 0) {
|
---|
| 210 | // be careful of using face (for location hinting) only when it exists
|
---|
| 211 | _supervertex[index].vertex = _TR.insert(Point(points_to_add[ia].first,
|
---|
| 212 | points_to_add[ia].second),face);}
|
---|
| 213 | else {
|
---|
| 214 | _supervertex[index].vertex = _TR.insert(Point(points_to_add[ia].first,
|
---|
| 215 | points_to_add[ia].second));
|
---|
| 216 | }
|
---|
| 217 | // we are not up to dealing with coincident vertices, so make
|
---|
| 218 | // sure the user knows!
|
---|
| 219 | _CrashIfVertexPresent(_supervertex[index].vertex, index);
|
---|
| 220 | _supervertex[index].vertex->info() = index;
|
---|
| 221 |
|
---|
| 222 | // first find nearest neighbour of "newpoint" (shorthand for
|
---|
| 223 | // _supervertex[index].vertex); while we're at it, for each of the
|
---|
| 224 | // voronoi neighbours, "D", of newpoint, examine whether newpoint is
|
---|
| 225 | // closer to "D" than D's current nearest neighbour -- when this
|
---|
| 226 | // occurs, put D into indices_of_updated_neighbours.
|
---|
| 227 | //
|
---|
| 228 | // manually put newpoint on indices_of_updated_neighbours
|
---|
| 229 | indices_of_updated_neighbours.push_back(index);
|
---|
| 230 | _SetAndUpdateNearest(index, indices_of_updated_neighbours);
|
---|
| 231 | }
|
---|
| 232 |
|
---|
| 233 | // for Voronoi neighbours j of any of the removed points for which
|
---|
| 234 | // one of those removed points was the nearest neighbour,
|
---|
| 235 | // redetermine the nearest neighbour of j and add j onto the vector
|
---|
| 236 | // of indices_of_updated_neighbours.
|
---|
| 237 | set<int>::iterator it2 = NeighbourUnion.begin();
|
---|
| 238 | for ( ; it2 != NeighbourUnion.end(); ++it2 ) {
|
---|
| 239 | int j = *it2;
|
---|
| 240 | // the if avoids the vertex at infinity, which gets a negative index
|
---|
| 241 | if( j != INFINITE_VERTEX ) {
|
---|
| 242 | // this is where we check if the nearest neighbour of j was one
|
---|
| 243 | // of the removed points
|
---|
| 244 | if (indices_removed.count(_supervertex[j].NNindex)) {
|
---|
| 245 | if (_verbose) cout << "j " << j << endl;
|
---|
| 246 | _SetNearest(j);
|
---|
| 247 | indices_of_updated_neighbours.push_back(j);
|
---|
| 248 | if (_verbose) cout << "NN of " << j << " : "
|
---|
| 249 | << _supervertex[j].NNindex
|
---|
| 250 | << ", dist = " << _supervertex[j].NNdistance <<endl;
|
---|
| 251 | }
|
---|
| 252 | }
|
---|
| 253 | }
|
---|
| 254 |
|
---|
| 255 | }
|
---|
| 256 |
|
---|
| 257 |
|
---|
| 258 | //----------------------------------------------------------------------
|
---|
| 259 | /// Determines the index and distance of the nearest neighbour to
|
---|
| 260 | /// point j and puts the information into the _supervertex entry for j.
|
---|
| 261 | void DnnPlane::_SetNearest (const int & j) {
|
---|
| 262 | Vertex_handle current = _supervertex[j].vertex;
|
---|
| 263 | Vertex_circulator vc = _TR.incident_vertices(current);
|
---|
| 264 | Vertex_circulator done = vc;
|
---|
| 265 | double dist;
|
---|
| 266 | double mindist = HUGE_DOUBLE; // change this to "HUGE" or max_double?
|
---|
| 267 | Vertex_handle nearest = _TR.infinite_vertex();
|
---|
| 268 |
|
---|
| 269 | // when there is only one finite point left in the triangulation,
|
---|
| 270 | // there are no triangles. Presumably this is why voronoi returns
|
---|
| 271 | // NULL for the incident vertex circulator. Check if this is
|
---|
| 272 | // happening before circulating over it... (Otherwise it crashes
|
---|
| 273 | // when looking for neighbours of last point)
|
---|
| 274 | if (vc != NULL) do {
|
---|
| 275 | if ( vc->info().val() != INFINITE_VERTEX) {
|
---|
| 276 | // find distance between j and its Voronoi neighbour (vc)
|
---|
| 277 | if (_verbose) cout << current->info().val() << " " << vc->info().val() << endl;
|
---|
| 278 | dist = _euclid_distance(current->point(), vc->point());
|
---|
| 279 | // check if j is closer to vc than vc's currently registered
|
---|
| 280 | // nearest neighbour (and update things if it is)
|
---|
| 281 | if (dist < mindist) {
|
---|
| 282 | mindist = dist; nearest = vc;
|
---|
| 283 | if (_verbose) cout << "nearer ";
|
---|
| 284 | }
|
---|
| 285 | if (_verbose) cout << vc->point() << "; "<< dist << endl;
|
---|
| 286 | }
|
---|
| 287 | } while (++vc != done); // move on to next Voronoi neighbour
|
---|
| 288 | // set j's supervertex info about nearest neighbour
|
---|
| 289 | _supervertex[j].NNindex = nearest->info().val();
|
---|
| 290 | _supervertex[j].NNdistance = mindist;
|
---|
| 291 | }
|
---|
| 292 |
|
---|
| 293 | //----------------------------------------------------------------------
|
---|
| 294 | /// Determines and stores the nearest neighbour of j, and where
|
---|
| 295 | /// necessary update the nearest-neighbour info of Voronoi neighbours
|
---|
| 296 | /// of j;
|
---|
| 297 | ///
|
---|
| 298 | /// For each voronoi neighbour D of j if the distance between j and D
|
---|
| 299 | /// is less than D's own nearest neighbour, then update the
|
---|
| 300 | /// nearest-neighbour info in D; push D's index onto
|
---|
| 301 | /// indices_of_updated_neighbours
|
---|
| 302 | ///
|
---|
| 303 | /// Note that j is NOT pushed onto indices_of_updated_neighbours --
|
---|
| 304 | /// if you want it there, put it there yourself.
|
---|
| 305 | ///
|
---|
| 306 | /// NB: note that we have _SetAndUpdateNearest as a completely
|
---|
| 307 | /// separate routine from _SetNearest because we want to
|
---|
| 308 | /// use one single ciruclation over voronoi neighbours to find the
|
---|
| 309 | /// nearest neighbour and to update the voronoi neighbours if need
|
---|
| 310 | /// be.
|
---|
| 311 | void DnnPlane::_SetAndUpdateNearest(
|
---|
| 312 | const int & j,
|
---|
| 313 | vector<int> & indices_of_updated_neighbours) {
|
---|
| 314 |
|
---|
| 315 | Vertex_handle current = _supervertex[j].vertex;
|
---|
| 316 | Vertex_circulator vc = _TR.incident_vertices(current);
|
---|
| 317 | Vertex_circulator done = vc;
|
---|
| 318 | double dist;
|
---|
| 319 | double mindist = HUGE_DOUBLE; // change this to "HUGE" or max_double?
|
---|
| 320 | Vertex_handle nearest = _TR.infinite_vertex();
|
---|
| 321 |
|
---|
| 322 | // when there is only one finite point left in the triangulation,
|
---|
| 323 | // there are no triangles. Presumably this is why voronoi returns
|
---|
| 324 | // NULL for the incident vertex circulator. Check if this is
|
---|
| 325 | // happening before circulating over it... (Otherwise it crashes
|
---|
| 326 | // when looking for neighbours of last point)
|
---|
| 327 | if (vc != NULL) do {
|
---|
| 328 | if (vc->info().val() != INFINITE_VERTEX) {
|
---|
| 329 | if (_verbose) cout << current->info().val() << " " << vc->info().val() << endl;
|
---|
| 330 | // find distance between j and its Voronoi neighbour (vc)
|
---|
| 331 | dist = _euclid_distance(current->point(), vc->point());
|
---|
| 332 | // update the mindist if we are closer than anything found so far
|
---|
| 333 | if (dist < mindist) {
|
---|
| 334 | mindist = dist; nearest = vc;
|
---|
| 335 | if (_verbose) cout << "nearer ";
|
---|
| 336 | }
|
---|
| 337 | // find index corresponding to vc for easy manipulation
|
---|
| 338 | int vcindx = vc->info().val();
|
---|
| 339 | if (_verbose) cout << vc->point() << "; "<< dist << endl;
|
---|
| 340 | // check if j is closer to vc than vc's currently registered
|
---|
| 341 | // nearest neighbour (and update things if it is)
|
---|
| 342 | if (dist < _supervertex[vcindx].NNdistance) {
|
---|
| 343 | if (_verbose) cout << vcindx << "'s NN becomes " << current->info().val() << endl;
|
---|
| 344 | _supervertex[vcindx].NNdistance = dist;
|
---|
| 345 | _supervertex[vcindx].NNindex = j;
|
---|
| 346 | indices_of_updated_neighbours.push_back(vcindx);
|
---|
| 347 | }
|
---|
| 348 | }
|
---|
| 349 | } while (++vc != done); // move on to next Voronoi neighbour
|
---|
| 350 | // set j's supervertex info about nearest neighbour
|
---|
| 351 | _supervertex[j].NNindex = nearest->info().val();
|
---|
| 352 | _supervertex[j].NNdistance = mindist;
|
---|
| 353 | }
|
---|
| 354 |
|
---|
| 355 | FASTJET_END_NAMESPACE
|
---|
| 356 |
|
---|
| 357 | #endif // DROP_CGAL
|
---|