Skip to Content.
Sympa Menu

cgal-discuss - Re: [cgal-discuss] Add information in Triangulation while inserting a range of points

Subject: CGAL users discussion list

List archive

Re: [cgal-discuss] Add information in Triangulation while inserting a range of points


Chronological Thread 
  • From: Marc Alexa <>
  • To:
  • Subject: Re: [cgal-discuss] Add information in Triangulation while inserting a range of points
  • Date: Wed, 15 May 2019 11:09:29 +0200
  • Authentication-results: mail2-smtp-roc.national.inria.fr; spf=None ; spf=Pass ; spf=None
  • Ironport-phdr: 9a23:kuDZ1xG2ImGB/IYExsFFoZ1GYnF86YWxBRYc798ds5kLTJ7zocmwAkXT6L1XgUPTWs2DsrQY0rOQ4vCrADRIyK3CmUhKSIZLWR4BhJdetC0bK+nBN3fGKuX3ZTcxBsVIWQwt1Xi6NU9IBJS2PAWK8TW94jEIBxrwKxd+KPjrFY7OlcS30P2594HObwlSizexfK5+IA+yoAjSucUanJduIbstxxXUpXdFZ/5Yzn5yK1KJmBb86Maw/Jp9/ClVpvks6c1OX7jkcqohVbBXAygoPG4z5M3wqBnMVhCP6WcGUmUXiRVHHQ7I5wznU5jrsyv6su192DSGPcDzULs5Vyiu47ttRRT1jioMKjw3/3zNisFojKxUvB2uqQFxzY7afo+aNvlwcKTGcNwAWWZBW9xcVyxdDo6+aYYEEuoPPfxfr4n4v1YAqgGxBROwC+jy1jJIgmH53KIg3O88FgzG2RYvH8gSv3jOttr1MLkdUO+vw6TTwjXDaulZ2Tb56ITSbh8hpvSMUKt2fMHMx0cvEAbFgU+RqYzjJz6V0P4CvHOA4OpkS+2jkXIoqwZ0ojS3x8csjJPJhoMPxVze+yV52p45KsG3SEFhZd6oCpxQtzuVN4ZwX8gsQHlotT4kxrEavZO3ZisHxZQ9yxLCdfCKcJKE7x3hWeqJIjp1i2hpdK+9ihqu60Ss1OPxW8iu3FtLsCZIlMTHuGoX2BzJ8MeHT+Nw/ke/1jaL0ADe8uRELlo1larfMpIgzLswmocKvUTNESL7ml/6jKCRdkUj9eio7/robq/6qZ+bMo94kgD+MqIwlcyjGek0LBQCUmyB9em/1LDv51D1TKtJg/EsnaTUsojWJcEBqa64Bw9V3Jwj6xG6Dzq+1dQXh2MHI05fdB2di4jmJV7PL+rjAPewhlSjijZrx/TcMrL9BZXNK2DPkK39crZl905c1A0zwMhD6JJbEL4BJOv/VVLwtNzDEhA5Lhe0w/38BdVm1oIeXHqPDbWDPKPTt1+I/OMvLPOWaI8bojauY8Uj/OPk2H8lhUcGL+7uxooScHn+H/J8Ikzfb2CrmcYECW5NvwwwS6vhh1SGFDJSfH2vRLlv2jZuA42vCcLPR5umnaea9Ca9BJxfIG5cWX6WFnK9UoyeUL8lbC+CK4c1lzUeXv6oT4Ix3DmhsQb7z/xsKe+CqX5Qjo7qyNUgv76brho17zEhV53AgVHIdHl9myYzfxFz3K17phYgmFKK0Kw9nPkBUNIPuLVGVQA1MZOaxOt/WYirC1DxO+yRQVPjee2IRDQ4T9Y/2dgLOh8vFNCrjxSF1C2vUeZMy+67Qacs+6eZ5EDfYt5nwi+fhqYkhlgiBMBIMD/+iw==

I think the question is less why one would like to call some function, but consistency among derived classes (see below). If you absolutely need some reasons:
1. One could hope that a triangulation may be generated slightly faster than Delaunay if less flips are performed, such as for the pushing/placing triangulation. 
(2. [much less of a reason] an order dependent triangulation would allow creating many different triangulations in a convenient way, just by permuting the order and then calling the constructor. I would have liked to make use of this.)

Consistency: Delaunay_triangulation_3 has the following methods:
std::ptrdiff_t CGAL::Delaunay_triangulation_3< Traits, TDS, LP, SLDS >::insert (PointInputIterator first, 
PointInputIterator last 
)
std::ptrdiff_t CGAL::Delaunay_triangulation_3< Traits, TDS, LP, SLDS >::insert (PointWithInfoInputIterator first, 
PointWithInfoInputIterator last 
)
The second one apparently accepts std::pair as a way to provide points with info. 

Triangulation_3, from which Delaunay_triangulation_3 is derived, offers
std::ptrdiff_t CGAL::Triangulation_3< Traits, TDS, SLDS >::insert (InputIterator first, 
InputIterator last 
)
and the InputIterator seems to only accept Points, but not pairs.  

On the other hand, the vertex type used for Delaunay_triangulation_3 is Triangulation_vertex_base_3 or Triangulation_vertex_base_with_info_3 (and not something specialized for Delaunay_triangulation_3). 

From the perspective of using these interfaces this may seem odd. 

BTW: what is even more inconvenient is that Regular_triangulation_3 has no _with_info type at all. 

It would really be nice to have an easy way to access the vertices in a triangulation based on some existing integer indices. One way, as Mael pointed out, is to insert the vertices one by one. This allows creating maps from and to the vertex handles. But it is inefficient. It seems the way out is to apply spatial sorting to the points, keep track of the permutation, and then insert one by one. Is this expected to be as fast as calling the insert function with a range (modulo parallelism)?

-Marc





Let me correct myself, The constructor with a range, as well as the insert function   
with a range exist, and they perform a spatial sort, so we do what I  explained to
better not do in my previous mail.   And I didn't adress at all your request for
dealing with points with info.

We can look into this, but can you tell me why you would construct a Triangulation_3?

best,

andreas

On 5/13/2019 10:10 AM, Andreas Fabri wrote:


We could add an insertion of a range to the Triangulation_3 class, but then no
reshuffling of the points should happen as the insertion order implies which
triangles you have.  This is different for the Delaunay/Regular triangulation
where the resulting triangulation is unique, that is independent from the
insertion order.

andreas


On 5/13/2019 9:01 AM, Marc Alexa wrote:
Yes, CGAL is sometimes mysterious in what types are accepted by what functions. I agree you’d expect that the insert method for Triangulation_3 should be able to understand how to handle the pairs if Delaunay_triangulation_3 does. But, no. If you insist on using the pairs and are ok with initially getting a Delaunay triangulation, the following code works just fine:

#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Triangulation_3.h>
#include <CGAL/Delaunay_Triangulation_3.h>
#include <CGAL/Triangulation_cell_base_3.h>
#include <CGAL/Triangulation_vertex_base_with_info_3.h>
#include <vector>
typedef CGAL::Exact_predicates_inexact_constructions_kernel         K; typedef CGAL::Triangulation_vertex_base_with_info_3<unsigned, K>    Vb; typedef CGAL::Triangulation_cell_base_3<K>                 Cb; typedef CGAL::Triangulation_data_structure_3<Vb, Cb>                Tds; //Use the Fast_location tag. Default or Compact_location works too. typedef CGAL::Triangulation_3<K, Tds> Triangulation; typedef CGAL::Delaunay_triangulation_3<K, Tds> DT; typedef Triangulation::Point Point; int main() { std::vector< std::pair<Point, unsigned> > points; points.push_back(std::make_pair(Point(0, 0, 0), 0)); points.push_back(std::make_pair(Point(1, 0, 0), 1)); points.push_back(std::make_pair(Point(0, 1, 0), 2)); points.push_back(std::make_pair(Point(0, 0, 1), 3)); points.push_back(std::make_pair(Point(2, 2, 2), 4)); points.push_back(std::make_pair(Point(-1, 0, 1), 5));
    Triangulation T = DT(points.begin(), points.end());
    CGAL_assertion(T.number_of_vertices() == 6);

return 0; }
I believe speed-wise there is very little to be gained from not using Delaunay. If you don’t care about what coordinate in space carries which number, you could insert the range of points first and then set the info later. In this case, however, the order of the points you get form the vertex iterator will likely differ from the order in which you provided the points. 


-Marc


On 13. May 2019, at 04:25, sanlingtonCGAL <> wrote:

Hi,

I know the Delaunay and Regular triangulation in CGAL can already add
information to vertices while inserting /a range of points/, as introduced
in
https://doc.cgal.org/latest/Triangulation_3/index.html#Triangulation_3Delaunay

But it seems that the Triangulation_2/3 (not Delaunay triangulation or
Regular triangulation) can't do that, am I correct?

What I want to do is to construct Triangulation by inserting a range of
points, and each point or vertex would carry extra information, like an
index.
I am using codes like:

#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Triangulation_3.h>
#include <CGAL/Triangulation_cell_base_3.h>
#include <CGAL/Triangulation_vertex_base_with_info_3.h>
#include <vector>
typedef CGAL::Exact_predicates_inexact_constructions_kernel         K;
typedef CGAL::Triangulation_vertex_base_with_info_3<unsigned, K>    Vb;
typedef CGAL::Triangulation_cell_base_3<K>                 Cb;
typedef CGAL::Triangulation_data_structure_3<Vb, Cb>                Tds;
//Use the Fast_location tag. Default or Compact_location works too.
typedef CGAL::Triangulation_3<K, Tds> Triangulation;
typedef Triangulation::Point                                            
Point;
int main()
{
std::vector< std::pair<Point, unsigned> > points;
points.push_back(std::make_pair(Point(0, 0, 0), 0));
points.push_back(std::make_pair(Point(1, 0, 0), 1));
points.push_back(std::make_pair(Point(0, 1, 0), 2));
points.push_back(std::make_pair(Point(0, 0, 1), 3));
points.push_back(std::make_pair(Point(2, 2, 2), 4));
points.push_back(std::make_pair(Point(-1, 0, 1), 5));

Triangulation T(points.begin(), points.end());
CGAL_assertion(T.number_of_vertices() == 6);

return 0;
}

The above codes are similar to those in the manual document, I just replaced
the Delaunay with Triangulation. But the compiler gave me errors like:

c:\program files (x86)\microsoft visual
studio\2017\community\vc\tools\msvc\14.16.27023\include\xmemory0(881): error
C2664: 'CGAL::Point_3<Kernel_>::Point_3(CGAL::Point_3<Kernel_> &&)': cannot
convert argument 1 from '_Ty' to 'const CGAL::Origin &'
1>        with
1>        [
1>            Kernel_=CGAL::Epick
1>        ]
1>        and
1>        [
1>            _Ty=std::pair<Point,unsigned int>
1>        ]
1>c:\program files (x86)\microsoft visual
studio\2017\community\vc\tools\msvc\14.16.27023\include\xmemory0(879): note:
Reason: cannot convert from '_Ty' to 'const CGAL::Origin'
1>        with
1>        [
1>            _Ty=std::pair<Point,unsigned int>
1>        ]


What should I do to add information to point/vertex of a Triangulation while
inserting a range of points?

Thank you.





--
Sent from: http://cgal-discuss.949826.n4.nabble.com/

--
You are currently subscribed to cgal-discuss.
To unsubscribe or access the archives, go to
https://sympa.inria.fr/sympa/info/cgal-discuss



-- 
Andreas Fabri, PhD
Chief Officer, GeometryFactory
Editor, The CGAL Project

phone: +33.492.954.912    skype: andreas.fabri
-- 
Andreas Fabri, PhD
Chief Officer, GeometryFactory
Editor, The CGAL Project

phone: +33.492.954.912    skype: andreas.fabri




Archive powered by MHonArc 2.6.18.

Top of Page