Fork me on GitHub

source: git/external/fastjet/DnnPlane.cc@ 45094d9

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 45094d9 was 1d208a2, checked in by Pavel Demin <pavel.demin@…>, 8 years ago

update FastJet library to 3.2.1 and Nsubjettiness library to 2.2.4

  • Property mode set to 100644
File size: 24.5 KB
Line 
1//FJSTARTHEADER
2// $Id: DnnPlane.cc 3917 2015-07-03 14:07:50Z salam $
3//
4// Copyright (c) 2005-2014, 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#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
38using namespace std;
39
40FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh
41
42const double DnnPlane::DISTANCE_FOR_CGAL_CHECKS=1.0e-12;
43
44
45/// Initialiser from a set of points on an Eta-Phi plane, where both
46/// eta and phi can have arbitrary ranges
47DnnPlane::DnnPlane(const vector<EtaPhi> & input_points,
48 const bool & verbose ) {
49
50 _verbose = verbose;
51 int n = input_points.size();
52
53 // construct Voronoi diagram in such a way as to get the vertex handles
54 // and remember to set CGAL info with the index of the vertex
55 SuperVertex sv;
56 for (int i = 0; i < n; i++) {
57 sv.vertex =
58 _TR.insert(Point(input_points[i].first, input_points[i].second));
59
60 // check if we are dealing with coincident vertices
61 int coinciding_index = _CheckIfVertexPresent(sv.vertex, i);
62 if (coinciding_index == i){
63 // we need to associate an index to each vertex -- thus when we get
64 // a vertex (e.g. as a nearest neighbour) from CGAL, we will be
65 // able to figure out which particle it corresponded to.
66 sv.vertex->info() = sv.coincidence = i;
67 } else {
68 //cout << " coincident with " << coinciding_index << endl;
69 // the new vertex points to the already existing one and we
70 // record the coincidence
71 //
72 // Note that we must not only set the coincidence of the
73 // currently-added particle, the one it coincides with also
74 // needs be updated (taking into account that it might already
75 // coincide with another one)
76 //
77 // An example may help. Say coinciding_index = i1 and we're adding i2==i.
78 // Then _sv[i2].coincidence = i1; _sv[i1].coincidence = i2. In both
79 // cases sv.vertex->info() == i1;
80 //
81 // Later on we add i3; we find out that its coinciding index is i1;
82 // so we set _sv[i3].coincidence = i2 and sv[i1].coincidence = i3.
83 //
84 // This gives us the structure
85 // _supervertex[i1].coincidence == in
86 // _supervertex[i2].coincidence == i1
87 // ...
88 // _supervertex[in].coincidence == in-1
89 //
90 sv.coincidence = _supervertex[coinciding_index].coincidence; // handles cases with previous coincidences
91 _supervertex[coinciding_index].coincidence = i;
92 }
93
94 _supervertex.push_back(sv);
95 }
96
97 // label infinite vertex info with negative index
98 _TR.infinite_vertex()->info() = INFINITE_VERTEX;
99
100 // set up the structure that holds nearest distances and neighbours
101 for (int j = 0; j < n; j++) {_SetNearest(j);}
102
103}
104
105
106//----------------------------------------------------------------------
107/// Crashes if the given vertex handle already exists. Otherwise
108/// it does the bookkeeping for future such tests
109int DnnPlane::_CheckIfVertexPresent(
110 const Vertex_handle & vertex, const int its_index) {
111 // vertices that do not have the same geometric position as any
112 // other vertex so far added have info().val() == NEW_VERTEX -- this
113 // is ensured by the InitialisedInt class, which forms the "info"
114 // part of our
115 // CGAL::Triangulation_vertex_base_with_info_2<InitialisedInt,K>
116 //
117 // If the vertex coincides with one that already exists, then
118 // info().val() it's info().val() will have been updated (in
119 // DNN:DNN) to be equal to a vertex "index".
120 if (vertex->info().val() != NEW_VERTEX) {
121 if (_crash_on_coincidence){
122 ostringstream err;
123 err << "Error: DnnPlane::_CheckIfVertexPresent"
124 << "Point "<<its_index<<" coincides with point "
125 <<vertex->info().val() << endl;
126 throw DnnError(err.str());
127 }
128 return vertex->info().val();
129 }
130
131 return its_index;
132}
133
134
135//----------------------------------------------------------------------
136/// remove the points labelled by the vector indices_to_remove, and
137/// add the points specified by the vector points_to_add
138/// (corresponding indices will be calculated automatically); the
139/// idea behind this routine is that the points to be added will
140/// somehow be close to the one or other of the points being removed
141/// and this can be used by the implementation to provide hints for
142/// inserting the new points in whatever structure it is using. In a
143/// kt-algorithm the points being added will be a result of a
144/// combination of the points to be removed -- hence the proximity
145/// is (more or less) guaranteed.
146void DnnPlane::RemoveAndAddPoints(
147 const vector<int> & indices_to_remove,
148 const vector<EtaPhi> & points_to_add,
149 vector<int> & indices_added,
150 vector<int> & indices_of_updated_neighbours) {
151
152 if (_verbose) cout << "Starting DnnPlane::RemoveAndAddPoints" << endl;
153
154 // build set of UNION of Voronoi neighbours of a pair of nearest
155 // neighbours
156 set<int> NeighbourUnion;
157 // later on it will be convenient to have access to a set (rather
158 // than vector) of indices being removed
159 set<int> indices_removed;
160
161 // for each of the indices to be removed add the voronoi
162 // neighbourhood to the NeighbourUnion set as well as the coinciding
163 // points that had the current point as coincidence before.
164 for (size_t ir = 0; ir < indices_to_remove.size(); ir++) {
165 int index = indices_to_remove[ir];
166 indices_removed.insert(index);
167 if (_verbose) cout << " scheduling point " << index << " for removal" << endl;
168
169 if (_supervertex[index].coincidence != index){
170 // we have a coincidence
171 //
172 // The only one of the coincident points that has to be
173 // inserted in the neighbourhood list (and thus updated) is the
174 // one that has 'index' as coincidence.
175 int new_index = _supervertex[index].coincidence;
176 while (_supervertex[new_index].coincidence != index)
177 new_index = _supervertex[new_index].coincidence;
178 if (_verbose) cout << " inserted coinciding " << new_index << " to neighbours union" << endl;
179 NeighbourUnion.insert(new_index);
180
181 // if this is the point among the coiciding ones that holds the
182 // CGAL vertex, then also insert the CGAL neighbours, otherwise
183 // just skip that step.
184 if (index != _supervertex[index].vertex->info().val()) continue;
185 }
186
187 // have a circulators that will go round the Voronoi neighbours of
188 // _supervertex[index1].vertex
189 Vertex_circulator vc = _TR.incident_vertices(_supervertex[index].vertex);
190 Vertex_circulator done = vc;
191 if (vc != NULL){ // a safety check in case there is no Voronoi
192 // neighbour (which may happen e.g. if we just
193 // have a bunch of coincident points)
194 do {
195 // if a neighbouring vertex is not the infinite vertex, then add it
196 // to our union of neighbouring vertices.
197 if (_verbose) cout << "examining " << vc->info().val() << endl;
198 if (vc->info().val() != INFINITE_VERTEX) {
199 // NB: from it=1 onwards occasionally it might already have
200 // been inserted -- but double insertion still leaves only one
201 // copy in the set, so there's no problem
202 NeighbourUnion.insert(vc->info().val());
203 if (_verbose) cout << " inserted " << vc->info().val() << " to neighbours union" << endl;
204 }
205 } while (++vc != done);
206 }
207 }
208
209 if (_verbose) {
210 set<int>::iterator it = NeighbourUnion.begin();
211 cout << "Union of neighbours of combined points" << endl;
212 for ( ; it != NeighbourUnion.end(); ++it ) {
213 cout << *it << endl;
214 }
215 }
216
217 // update set, triangulation and supervertex info
218 for (size_t ir = 0; ir < indices_to_remove.size(); ir++) {
219 int index = indices_to_remove[ir];
220 if (_verbose) cout << " removing " << index << endl;
221
222 // NeighbourUnion should not contain the points to be removed
223 // (because later we will assume they still exist).
224 NeighbourUnion.erase(indices_to_remove[ir]);
225
226 // first deal with coincidences
227 if (_supervertex[index].coincidence != index){
228 int new_index = _supervertex[index].coincidence;
229
230 // if this is the point among the coiciding ones that "owns" the
231 // CGAL vertex we need to re-label the CGAL vertex so that it
232 // points to the coincident particle and set the current one to
233 // NULL
234 //
235 // This can be done only on the first point as they all share
236 // the same value
237 //
238 // Note that this has to be done before the following step since
239 // it will alter the coincidence information
240 if (index == _supervertex[index].vertex->info().val())
241 _supervertex[new_index].vertex->info() = new_index;
242
243 // we need to browse the coincidences until we end the loop, at
244 // which point we reset the coincidence of the point that has
245 // the current one as a coincidence
246 while (_supervertex[new_index].coincidence != index)
247 new_index = _supervertex[new_index].coincidence;
248 _supervertex[new_index].coincidence = _supervertex[index].coincidence;
249
250 // remove the coincidence on the point being removed and mark it
251 // as removed
252 _supervertex[index].coincidence = index;
253 _supervertex[index].vertex = NULL;
254
255 continue;
256 }
257
258 // points to be removed should also be eliminated from the
259 // triangulation and the supervertex structure should be updated
260 // to reflect the fact that the points are no longer valid.
261 _TR.remove(_supervertex[index].vertex);
262 if (_verbose) cout << "DnnPlane about to set _supervertex["<< index<<"].vertex to NULL" << endl;
263 _supervertex[index].vertex = NULL;
264 if (_verbose) cout << " value is " << (_is_not_null(_supervertex[index].vertex)) << endl;
265 }
266
267 // add new point: give a "hint" to the inserter that
268 // the new point should be added close to old points -- the easiest way
269 // of getting this is to take a point from the NeighbourUnion AFTER we have
270 // removed point1, point2, and to get one of its incident faces.
271 //
272 // This hinting improves speed by c. 25% for 10^4 points because it
273 // avoids the need for a costly (sqrt{N}) location step (at least
274 // with a non-hierarchical triangulation -- with a hierarchical one,
275 // this step could be done away with, though there will still be a
276 // cost of O(ln N) to pay.
277 //
278 // For some reason inserting the point before the two removals
279 // slows things down by c. 25%. This importance of the order
280 // is not understood.
281 //
282 // At some point it might be worth trying to select the "nearest"
283 // of the various points in the neighbour union to avoid large
284 // steps in cases where we have 0..2pi periodicity and the first member
285 // of the neighbour union happens to be on the wrong side.
286 Face_handle face;
287 //if (indices_to_remove.size() > 0) { // GS: use NeighbourUnion instead
288 // (safe also in case of coincidences)
289 if (NeighbourUnion.size() > 0) {
290 // face can only be found if there were points to remove in first place
291 face = _TR.incident_faces(
292 _supervertex[*NeighbourUnion.begin()].vertex);}
293 // make sure the output arrays are empty
294 indices_added.clear();
295 indices_of_updated_neighbours.clear();
296 for (size_t ia = 0; ia < points_to_add.size(); ia++) {
297 SuperVertex sv;
298 _supervertex.push_back(sv);
299 int index = _supervertex.size()-1;
300 indices_added.push_back(index);
301 if (_verbose) cout << " adding " << index << " at "
302 << points_to_add[ia].first<< " " << points_to_add[ia].second << endl;
303
304 //if (indices_to_remove.size() > 0) {
305 if (NeighbourUnion.size() > 0) {
306 // be careful of using face (for location hinting) only when it exists
307 _supervertex[index].vertex = _TR.insert(Point(points_to_add[ia].first,
308 points_to_add[ia].second),face);}
309 else {
310 _supervertex[index].vertex = _TR.insert(Point(points_to_add[ia].first,
311 points_to_add[ia].second));
312 }
313
314 // check if this leads to a coincidence
315 int coinciding_index = _CheckIfVertexPresent(_supervertex[index].vertex, index);
316 if (coinciding_index == index){
317 // we need to associate an index to each vertex -- thus when we get
318 // a vertex (e.g. as a nearest neighbour) from CGAL, we will be
319 // able to figure out which particle it corresponded to.
320 _supervertex[index].vertex->info() = _supervertex[index].coincidence = index;
321 } else {
322 if (_verbose) cout << " coinciding with vertex " << coinciding_index << endl;
323 // the new vertex points to an already existing one and we
324 // record the coincidence
325 //
326 // we also update the NN of the coinciding particle (to avoid
327 // having to loop over the list of coinciding neighbours later)
328 // This is done first as it allows us to check if this is a new
329 // coincidence or a coincidence added to a particle that was
330 // previously "alone"
331 _supervertex[coinciding_index].NNindex = index;
332 _supervertex[coinciding_index].NNdistance = 0.0;
333 indices_of_updated_neighbours.push_back(coinciding_index);
334
335 // Note that we must not only set the coincidence of the
336 // currently-added particle, the one it coincides with also
337 // needs be updated (taking into account that it might already
338 // coincide with another one)
339 _supervertex[index].coincidence = _supervertex[coinciding_index].coincidence; // handles cases with previous coincidences
340 _supervertex[coinciding_index].coincidence = index;
341
342 }
343
344 // first find nearest neighbour of "newpoint" (shorthand for
345 // _supervertex[index].vertex); while we're at it, for each of the
346 // voronoi neighbours, "D", of newpoint, examine whether newpoint is
347 // closer to "D" than D's current nearest neighbour -- when this
348 // occurs, put D into indices_of_updated_neighbours.
349 //
350 // manually put newpoint on indices_of_updated_neighbours
351 indices_of_updated_neighbours.push_back(index);
352 _SetAndUpdateNearest(index, indices_of_updated_neighbours);
353
354 //cout << "Added: " << setprecision(20) << " ("
355 // << points_to_add[ia].first << "," << points_to_add[ia].second
356 // << ") with index " << index << endl;
357 }
358
359 // for Voronoi neighbours j of any of the removed points for which
360 // one of those removed points was the nearest neighbour,
361 // redetermine the nearest neighbour of j and add j onto the vector
362 // of indices_of_updated_neighbours.
363 set<int>::iterator it2 = NeighbourUnion.begin();
364 for ( ; it2 != NeighbourUnion.end(); ++it2 ) {
365 int j = *it2;
366 // the if avoids the vertex at infinity, which gets a negative index
367 if( j != INFINITE_VERTEX ) {
368 // this is where we check if the nearest neighbour of j was one
369 // of the removed points
370 if (indices_removed.count(_supervertex[j].NNindex)) {
371 if (_verbose) cout << "j " << j << endl;
372 _SetNearest(j);
373 indices_of_updated_neighbours.push_back(j);
374 if (_verbose) cout << "NN of " << j << " : "
375 << _supervertex[j].NNindex
376 << ", dist = " << _supervertex[j].NNdistance <<endl;
377 }
378 }
379 }
380
381 if (_verbose) cout << "Leaving DnnPlane::RemoveAndAddPoints" << endl;
382}
383
384//----------------------------------------------------------------------
385/// Determines the index and distance of the nearest neighbour to
386/// point j and puts the information into the _supervertex entry for j.
387void DnnPlane::_SetNearest (const int j) {
388 // first deal with the cases where we have a coincidence
389 if (_supervertex[j].coincidence != j){
390 _supervertex[j].NNindex = _supervertex[j].coincidence;
391 _supervertex[j].NNdistance = 0.0;
392 return;
393 }
394
395 // The code below entirely uses CGAL distance comparisons to compute
396 // the nearest neighbour. It has the mais drawback to induice a
397 // 10-20% time penalty so we switched to our own comparison (which
398 // only turns to CGAL for dangerous situations)
399 //
400 // Vertex_handle current = _supervertex[j].vertex;
401 // Vertex_circulator vc = _TR.incident_vertices(current);
402 // Vertex_circulator done = vc;
403 // Vertex_handle nearest = _TR.infinite_vertex();
404 // double mindist = HUGE_DOUBLE;
405 //
406 // // when there is only one finite point left in the triangulation,
407 // // there are no triangles. Presumably this is why voronoi returns
408 // // NULL for the incident vertex circulator. Check if this is
409 // // happening before circulating over it... (Otherwise it crashes
410 // // when looking for neighbours of last point)
411 // if (vc != NULL){
412 // // initialise the nearest vertex handle to the first incident
413 // // vertex that is not INFINITE_VERTEX
414 // while (vc->info().val() == INFINITE_VERTEX){
415 // vc++;
416 // if (vc==done) break; // if vc==done, then INFINITE_VERTEX is the
417 // // only element in the neighbourhood
418 // }
419 //
420 // // if there is just the infinite vertex, we have vc->info().val()
421 // // == INFINITE_VERTEX and nothing has to be done
422 // // otherwise, use the current vc as an initialisation
423 // if (vc->info().val() != INFINITE_VERTEX){
424 // nearest = vc; // initialisation to the first non-infinite vertex
425 //
426 // // and loop over the following ones
427 // while (++vc != done){
428 // // we should not compare with the infinite vertex
429 // if (vc->info().val() == INFINITE_VERTEX) continue;
430 //
431 // if (_verbose) cout << current->info().val() << " " << vc->info().val() << endl;
432 // // use CGAL's distance comparison to check if 'vc' is closer to
433 // // 'current' than the nearest so far (we include the == case for
434 // // safety though it should not matter in this precise case)
435 // if (CGAL::compare_distance_to_point(current->point(), vc->point(), nearest->point())!=CGAL::LARGER){
436 // nearest = vc;
437 // if (_verbose) cout << "nearer";
438 // }
439 // }
440 //
441 // // now compute the distance
442 // //
443 // // Note that since we're always using CGAL to compare distances
444 // // (and never the distance computed using _euclid_distance) we
445 // // should not worry about rounding errors in mindist
446 // mindist = _euclid_distance(current->point(), nearest->point());
447 // }
448 // }
449 //
450 // // set j's supervertex info about nearest neighbour
451 // _supervertex[j].NNindex = nearest->info().val();
452 // _supervertex[j].NNdistance = mindist;
453
454 Vertex_handle current = _supervertex[j].vertex;
455 Vertex_circulator vc = _TR.incident_vertices(current);
456 Vertex_circulator done = vc;
457 double dist;
458 double mindist = HUGE_DOUBLE; // change this to "HUGE" or max_double?
459 Vertex_handle nearest = _TR.infinite_vertex();
460
461 // when there is only one finite point left in the triangulation,
462 // there are no triangles. Presumably this is why voronoi returns
463 // NULL for the incident vertex circulator. Check if this is
464 // happening before circulating over it... (Otherwise it crashes
465 // when looking for neighbours of last point)
466 if (vc != NULL) do {
467 if ( vc->info().val() != INFINITE_VERTEX) {
468 // find distance between j and its Voronoi neighbour (vc)
469 if (_verbose) cout << current->info().val() << " " << vc->info().val() << endl;
470
471 // check if j is closer to vc than vc's currently registered
472 // nearest neighbour (and update things if it is)
473 if (_is_closer_to(current->point(), vc->point(), nearest, dist, mindist)){
474 nearest = vc;
475 if (_verbose) cout << "nearer ";
476 }
477 if (_verbose) cout << vc->point() << "; "<< dist << endl;
478 }
479 } while (++vc != done); // move on to next Voronoi neighbour
480
481 // set j's supervertex info about nearest neighbour
482 _supervertex[j].NNindex = nearest->info().val();
483 _supervertex[j].NNdistance = mindist;
484}
485
486//----------------------------------------------------------------------
487/// Determines and stores the nearest neighbour of j, and where
488/// necessary updates the nearest-neighbour info of Voronoi neighbours
489/// of j;
490///
491/// For each voronoi neighbour D of j if the distance between j and D
492/// is less than D's own nearest neighbour, then update the
493/// nearest-neighbour info in D; push D's index onto
494/// indices_of_updated_neighbours
495///
496/// Note that j is NOT pushed onto indices_of_updated_neighbours --
497/// if you want it there, put it there yourself.
498///
499/// NB: note that we have _SetAndUpdateNearest as a completely
500/// separate routine from _SetNearest because we want to
501/// use one single circulation over voronoi neighbours to find the
502/// nearest neighbour and to update the voronoi neighbours if need
503/// be.
504void DnnPlane::_SetAndUpdateNearest(
505 const int j,
506 vector<int> & indices_of_updated_neighbours) {
507 //cout << "SetAndUpdateNearest for point " << j << endl;
508 // first deal with coincidences
509 if (_supervertex[j].coincidence != j){
510 _supervertex[j].NNindex = _supervertex[j].coincidence;
511 _supervertex[j].NNdistance = 0.0;
512 //cout << " set to coinciding point " << _supervertex[j].coincidence << endl;
513 return;
514 }
515
516 Vertex_handle current = _supervertex[j].vertex;
517 Vertex_circulator vc = _TR.incident_vertices(current);
518 Vertex_circulator done = vc;
519 double dist;
520 double mindist = HUGE_DOUBLE; // change this to "HUGE" or max_double?
521 Vertex_handle nearest = _TR.infinite_vertex();
522
523 // when there is only one finite point left in the triangulation,
524 // there are no triangles. Presumably this is why voronoi returns
525 // NULL for the incident vertex circulator. Check if this is
526 // happening before circulating over it... (Otherwise it crashes
527 // when looking for neighbours of last point)
528 if (vc != NULL) do {
529 if (vc->info().val() != INFINITE_VERTEX) {
530 if (_verbose) cout << current->info().val() << " " << vc->info().val() << endl;
531
532 // update the mindist if we are closer than anything found so far
533 if (_is_closer_to(current->point(), vc->point(), nearest, dist, mindist)){
534 nearest = vc;
535 if (_verbose) cout << "nearer ";
536 }
537
538 // find index corresponding to vc for easy manipulation
539 int vcindx = vc->info().val();
540 if (_verbose) cout << vc->point() << "; "<< dist << endl;
541
542 if (_is_closer_to_with_hint(vc->point(), current->point(),
543 _supervertex[_supervertex[vcindx].NNindex].vertex,
544 dist, _supervertex[vcindx].NNdistance)){
545 if (_verbose) cout << vcindx << "'s NN becomes " << current->info().val() << endl;
546 _supervertex[vcindx].NNindex = j;
547 indices_of_updated_neighbours.push_back(vcindx);
548 }
549
550 // original code without the use of CGAL distance in potentially
551 // dangerous cases
552 //
553 // // check if j is closer to vc than vc's currently registered
554 // // nearest neighbour (and update things if it is)
555 // //
556 // // GS: originally, the distance test below was a strict <. It
557 // // has to be <= because if the two distances are ==, it is
558 // // possible that the old NN is no longer connected to vc in
559 // // the triangulation, and we are sure that the newly
560 // // inserted point (j) is (since we loop over j's
561 // // neighbouring points in the triangulation).
562 // if (dist <= _supervertex[vcindx].NNdistance) {
563 // if (_verbose) cout << vcindx << "'s NN becomes " << current->info().val() << endl;
564 // _supervertex[vcindx].NNdistance = dist;
565 // _supervertex[vcindx].NNindex = j;
566 // indices_of_updated_neighbours.push_back(vcindx);
567 // }
568 }
569 } while (++vc != done); // move on to next Voronoi neighbour
570 // set j's supervertex info about nearest neighbour
571 //cout << " set to point " << nearest->info().val() << endl;
572 _supervertex[j].NNindex = nearest->info().val();
573 _supervertex[j].NNdistance = mindist;
574}
575
576FASTJET_END_NAMESPACE
577
578#endif // DROP_CGAL
Note: See TracBrowser for help on using the repository browser.