Skip to Content.
Sympa Menu

cgal-discuss - Re: [cgal-discuss] Using finger search with Delaunay Triangulation

Subject: CGAL users discussion list

List archive

Re: [cgal-discuss] Using finger search with Delaunay Triangulation


Chronological Thread 
  • From: Damian Sheehy <>
  • To:
  • Subject: Re: [cgal-discuss] Using finger search with Delaunay Triangulation
  • Date: Sun, 19 Apr 2009 07:54:25 -0400
  • Domainkey-signature: a=rsa-sha1; c=nofws; d=gmail.com; s=gamma; h=mime-version:in-reply-to:references:date:message-id:subject:from:to :content-type; b=CMoHzzJMWKmWfDm9OkZhkO2bl9h/VMQneQHo1F/6hv6az26/7yd1BEYvZ4y8NZQsmR 66Rf58GlG7pRDRCvpZr7ThYRsi6O0Lo7cJzeVROb3bvs5Ku0uswxEMfvA4h5y8SH54lA CO0ViNcw79ZNSKqM7u3KkA/gVcugThQL9Pq08=

Hi Manuel,
 
For a random set of points the Delaunay hierarchy should build a triangulation in much less time than the conventional Delaunay without a hierarchy. The numbers do not reflect that, so something is wrong. I will run your example tomorrow and see what I get. . .
Regards,
 
Damian
 
 
On Sun, Apr 19, 2009 at 7:08 AM, Manuel Holtgrewe <> wrote:
Hi Damian,

Sorry, this was me being stupid. I was accidently using a strange
patch (I did not write myself) to Triangulation_hierarchy_2.

However, I really appreciate your comprehensive answers. I also had a
look at Olivier Deviller's paper.

I am confident that I am now compiling against the official 3.4
release of CGAL. For some reason, however, the delaunay triangulation
is very close to the delaunay hierarchy in the queries and the
construction. I attached the benchmark program's source file for any
debunking.

Do the times make sense or would you expect the delaunay hierarchy
being faster than the delaunay triangulation?

Delaunay Hierarchy:  finger query vs. non-finger query
        n       dt_build   dt_no_finger      dt_finger     dt_speedup
     dh_build   dh_no_finger      dh_finger     dh_speedup
        2       0.000011       0.404763       0.403647       1.002765
     0.000002       0.405436       0.405328       1.000266
        4       0.000009       0.089347       0.088662       1.007728
     0.000004       0.088971       0.088810       1.001812
        8       0.000012       0.041551       0.040761       1.019384
     0.000007       0.041474       0.040800       1.016520
       16       0.000025       0.164201       0.158971       1.032899
     0.000016       0.165585       0.160771       1.029943
       32       0.000032       0.069689       0.058045       1.200605
     0.000031       0.069807       0.058301       1.197357
       64       0.000068       0.089519       0.105669       0.847164
     0.000070       0.089585       0.105138       0.852070
      128       0.000138       0.117457       0.084255       1.394065
     0.000136       0.092525       0.064825       1.427303
      256       0.000286       0.118400       0.064962       1.822601
     0.000290       0.118712       0.064920       1.828590
      512       0.000559       0.140097       0.064911       2.158288
     0.000572       0.140355       0.064893       2.162870
     1024       0.001152       0.184033       0.081685       2.252957
     0.001178       0.184366       0.081764       2.254856
     2048       0.002357       0.249846       0.067586       3.696715
     0.002420       0.249937       0.067444       3.705851
     4096       0.004706       0.339099       0.073331       4.624231
     0.004792       0.339403       0.073350       4.627173
     8192       0.009519       0.470756       0.073420       6.411840
     0.009998       0.470833       0.073425       6.412451
    16384       0.019473       0.639252       0.066725       9.580394
     0.020252       0.640028       0.066687       9.597512
    32768       0.040319       0.907107       0.070604      12.847796
     0.042462       0.899606       0.067399      13.347467
    65536       0.082370       1.240469       0.067532      18.368593
     0.087055       1.242859       0.069638      17.847420
   131072       0.169797       1.758881       0.070463      24.961786
     0.178655       1.756205       0.069502      25.268456
   262144       0.345318       2.509199       0.071552      35.068168
     0.362546       2.516872       0.072444      34.742214
   524288       0.715394       3.665781       0.071931      50.962377
     0.750623       3.674317       0.071780      51.188616
  1048576       1.466444       6.713272       0.072994      91.970195
     1.532415       6.730648       0.072908      92.317098
  2097152       3.033258      15.357660       0.075456     203.531586
     3.172991      15.327048       0.075358     203.390026
  4194304       6.165663      30.508426       0.078001     391.128543
     6.438287      30.560829       0.078401     389.802267
  8388608      12.615620      52.426141       0.081205     645.603242
    13.192613      52.019313       0.081088     641.516262


-- Manuel



2009/4/19 Damian Sheehy <>:
> Hi Manuel,
>
> The Delaunay hierarchy is very fast search structure in and of itself,
> particularly if a start hint cannot be provided for each query. The purpose
> of the search hierarchy is to compute a good start location from which to
> walk. In your case you already have a very good starting location
> (Face_handle) for each query point. So why would you expect the search
> hierarchy to do anything for you? (Take a look at the paper: Google Olivier
> Devillers Delaunay hierarchy).
>
> You are not really leveraging functionality of the search hierarchy, but you
> are paying the overhead of building and maintaining it. Try your search
> without the hierarchy, providing the start hint as before. If that's not
> fast enough consider exploiting the structure of your problem; you know that
> your query points are distributed along a straight line. CGAL provides a
> very powerful function called a line walk that allows you to compute the set
> of triangles intersected by a line segment. It may help to take a close look
> at that. You may be able to set up an interval list from the line walk and
> replace your N queries by a single line walk plus a subsequent search for an
> interval.
>
>
> All the best!
>
> Damian
>
> On Sat, Apr 18, 2009 at 12:16 PM, Manuel Holtgrewe <>
> wrote:
>>
>> Hi Damian,
>>
>> Thank you for your answer.
>>
>> I have written a little program with which I am trying to measure the
>> speedup that the fingers give me. The source code is attached, I
>> tested it with CGAL 3.4 on a machine with Intel Xeon processors
>> @2.3Ghz.
>>
>> I get the following output:
>>
>> Delaunay Hierarchy:  finger query vs. non-finger query
>>         n      build  no_finger     finger    speedup
>>         2   0.000012   0.424343   0.425491   0.997302
>>         4   0.000018   0.086763   0.086414   1.004036
>>         8   0.000015   0.100489   0.099627   1.008651
>>        16   0.000033   0.114695   0.105578   1.086352
>>        32   0.000066   0.128658   0.123606   1.040871
>>        64   0.000134   0.135320   0.116303   1.163512
>>       128   0.000270   0.178627   0.143877   1.241526
>>       256   0.000539   0.174100   0.122833   1.417371
>>       512   0.001070   0.168988   0.150091   1.125904
>>      1024   0.002091   0.173133   0.160160   1.080999
>>      2048   0.004395   0.185442   0.175380   1.057372
>>      4096   0.008549   0.207307   0.195271   1.061638
>>      8192   0.017413   0.216792   0.208223   1.041154
>>     16384   0.033093   0.246830   0.232705   1.060700
>>     32768   0.067720   0.198592   0.188638   1.052769
>>     65536   0.135914   0.207161   0.195284   1.060819
>>    131072   0.281401   0.243547   0.232596   1.047081
>>    262144   0.561701   0.254221   0.240915   1.055233
>>    524288   1.137475   0.238146   0.225836   1.054509
>>   1048576   2.287607   0.250391   0.238223   1.051078
>>   2097152   4.772575   0.252274   0.239942   1.051396
>>   4194304   9.628133   0.282821   0.270887   1.044055
>>   8388608  19.543273   0.305526   0.292194   1.045627
>>
>> As you can see, I could only measure very small speedups. I would
>> expect them to be much bigger. Am I doing something wrong?
>>
>> -- Manuel
>>
>>
>>
>> 2009/4/17 Damian Sheehy <>:
>> > Hi Manuel,
>> >
>> > Based on visual inspection, your code looks correct. You do not have
>> > to special case the first nearest_vertex search; just pass a default
>> > Face_handle in the first iteration of the loop and start the loop from
>> > zero.
>> >
>> > Regards,
>> >
>> > Damian
>> >
>> > On 4/17/09, Manuel Holtgrewe <> wrote:
>> >> Hi,
>> >>
>> >> I am trying to perform finger search with the delaunay hierarchy code
>> >> from CGAL. That is, I have a random point set which is added as the
>> >> points of a delaunay hierarchy / triangulation. Then, I have a list of
>> >> query points which lie along a line through this point set. I want to
>> >> query for the nearest neighbours of the points on the line.
>> >>
>> >> This sequence of queries should offer good locality, i.e. if I know
>> >> the nearest neighbour x of point i on the line, the nearest neighbour
>> >> y of point i+1 should be very close to x. Thus, I am trying to use the
>> >> second parameter to Delaunay_triangulation_2::nearest_vertex(const
>> >> Point& p, Face_handle f= Face_handle()) const.
>> >>
>> >> Will the following code do this, i.e. when calling x->face() (where x
>> >> is the nearest neighbour of the i-th query point i), will this return
>> >> a face adjacent to x and this helps to speed up the queries?
>> >>
>> >> Thanks,
>> >> Manuel
>> >>
>> >> // typedefs/aliases
>> >> struct K : CGAL::Exact_predicates_inexact_constructions_kernel {};
>> >>
>> >> typedef CGAL::Triangulation_vertex_base_2<K>              Vbb;
>> >> typedef CGAL::Triangulation_hierarchy_vertex_base_2<Vbb>  Vb;
>> >> typedef CGAL::Triangulation_face_base_2<K>                Fb;
>> >> typedef CGAL::Triangulation_data_structure_2<Vb,Fb>       Tds;
>> >>
>> >> typedef CGAL::Delaunay_triangulation_2<K, Tds>            Dt;
>> >> typedef CGAL::Triangulation_hierarchy_2<Dt>
>> >> Triangulation_hierarchy;
>> >> typedef Triangulation_hierarchy::Vertex_circulator
>> >> Vertex_circulator_hierarchy;
>> >> typedef Triangulation_hierarchy::Point
>> >>  Point_hierarchy;
>> >> typedef Triangulation_hierarchy::Vertex
>> >> Vertex_hierarchy;
>> >> typedef Triangulation_hierarchy::Vertex_handle
>> >> Vertex_hierarchy_handle;
>> >> typedef Triangulation_hierarchy::Face_handle
>> >> Face_hierarchy_handle;
>> >>
>> >> void foo(const std::vector<Point_hierarchy> &points, const
>> >> std::vector<Point_hierarchy> &queries) {
>> >>  Triangulation_hierarchy dh;
>> >>   dh.insert(points.begin(), points.end());
>> >>   Vertex_hierarchy_handle v = dh.nearest_vertex(queries[0]);
>> >>   Face_hierarchy_handle f = v->face();
>> >>   for (int j = 1, jend = queries.size(); j < jend; ++j) {
>> >>     v = dh.nearest_vertex(queries[j], f);
>> >>     f = v->face();
>> >>   }
>> >> }
>> >> --
>> >> You are currently subscribed to cgal-discuss.
>> >> To unsubscribe or access the archives, go to
>> >> https://lists-sop.inria.fr/wws/info/cgal-discuss
>> >>
>> > --
>> > You are currently subscribed to cgal-discuss.
>> > To unsubscribe or access the archives, go to
>> > https://lists-sop.inria.fr/wws/info/cgal-discuss
>> >
>> --
>> You are currently subscribed to cgal-discuss.
>> To unsubscribe or access the archives, go to
>> https://lists-sop.inria.fr/wws/info/cgal-discuss
>
>
--
You are currently subscribed to cgal-discuss.
To unsubscribe or access the archives, go to
https://lists-sop.inria.fr/wws/info/cgal-discuss




Archive powered by MHonArc 2.6.16.

Top of Page