Subject: CGAL users discussion list
List archive
RE: [cgal-discuss] Delaunay triangulation - find closest point within convex hull
Chronological Thread
- From: Ungermann, Jörn <>
- To: "" <>
- Subject: RE: [cgal-discuss] Delaunay triangulation - find closest point within convex hull
- Date: Wed, 7 Sep 2011 08:43:34 +0200
- Accept-language: en-US, de-DE
- Acceptlanguage: en-US, de-DE
Hi Bernd,
thanks, this might be what I am looking for!
After fighting with CGAL for a bit and using Gmpzf for the polytope distance
(otherwise it gave me a runtime error when calculating the realizing point),
I get indeed the same point I got before (as my geometry works out so that
the nearest vertex is indeed part of the closest edge).
I currently use the following (reduced) code:
<<< snippet start
typedef CGAL::Exact_predicates_inexact_constructions_kernel Kernel;
typedef CGAL::Triangulation_vertex_base_with_info_2<unsigned, Kernel> Vb;
typedef CGAL::Triangulation_data_structure_2<Vb> Tds;
typedef CGAL::Delaunay_triangulation_2<Kernel, Tds> Delaunay;
typedef CGAL::Polytope_distance_d<CGAL::Polytope_distance_d_traits_2<Kernel,
CGAL::Gmpzf, double> > PolyDist;
double x, y;
Delaunay T;
Kernel::Point_2 p[1] = {Kernel::Point_2(x, y)};
CGAL::Polytope_distance_d<CGAL::Polytope_distance_d_traits_2<Kernel,
CGAL::Gmpzf, double> > poly_dist(
p, p + 1, T.points_begin(), T.points_end());
Kernel::Point_2
p2(CGAL::to_double(*poly_dist.realizing_point_q_coordinates_begin ()) /
CGAL::to_double(*(poly_dist.realizing_point_q_coordinates_begin () + 2)),
CGAL::to_double(*(poly_dist.realizing_point_q_coordinates_begin () + 1)) /
CGAL::to_double(*(poly_dist.realizing_point_q_coordinates_begin () + 2)));
std::vector<std::pair<Delaunay::Vertex_handle, Kernel::FT> > vertices;
Kernel::FT norm2 = CGAL::natural_neighbor_coordinates_vertex_2(T, p2,
std::back_inserter(vertices)).second;
snippet end >>>
The point p2 looks good, but the natural_neighbour stuff still doesn't give
me vertices and weights.
Is this a problem with the exact/inexact coordinate stuff? Or shouldn't the
natural neighbor interpolation work for points on the border of the convex
hull?
Regards,
Joern
PS: How can I get the realizing_point without referring to the coordinates
iterator? Using the function for that wouldn't compile with some Gmpzf to
double conversion problem.
> -----Original Message-----
> From: Bernd Gaertner
> [mailto:]
> Sent: Dienstag, 6. September 2011 12:26
> To:
>
> Subject: Re: [cgal-discuss] Delaunay triangulation - find closest point
> within convex hull
>
> Ungermann wrote:
> > So, is there a function, which calculates for a given point the
> nearest
> > point from within the convex hull of a triangulation?
>
> There is Polytope_distance_d, see
> http://www.cgal.org/Manual/latest/doc_html/cgal_manual/Polytope_distanc
> e_d_ref/Class_Polytope_distance_d.html#Cross_link_anchor_1763.
>
> You need to feed it with the given point on the one hand and the convex
> hull vertices on the other hand. What you get out is the closest point,
> along with the convex hull points that generate it. Note that the
> closest point may not be on a face incident to the nearest vertex, so
> the heuristic you currently use may fail. If it happens to work in your
> case, it is sufficient to feed the nearest vertex plus its incident
> vertices to Polytope_distance_d.
>
> Best,
> Bernd.
>
> --
> You are currently subscribed to cgal-discuss.
> To unsubscribe or access the archives, go to
> https://lists-sop.inria.fr/wws/info/cgal-discuss
Attachment:
smime.p7s
Description: S/MIME cryptographic signature
- [cgal-discuss] Delaunay triangulation - find closest point within convex hull, Ungermann , Jörn, 09/06/2011
- Re: [cgal-discuss] Delaunay triangulation - find closest point within convex hull, Bernd Gaertner, 09/06/2011
- RE: [cgal-discuss] Delaunay triangulation - find closest point within convex hull, Ungermann , Jörn, 09/07/2011
- Re: [cgal-discuss] Delaunay triangulation - find closest point within convex hull, Andreas Fabri, 09/06/2011
- Re: [cgal-discuss] Delaunay triangulation - find closest point within convex hull, Bernd Gaertner, 09/06/2011
Archive powered by MHonArc 2.6.16.