1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203
|
// Copyright (c) 2005-2009 INRIA Sophia-Antipolis (France).
// All rights reserved.
//
// This file is part of CGAL (www.cgal.org)
//
// $URL$
// $Id$
// SPDX-License-Identifier: LGPL-3.0-or-later OR LicenseRef-Commercial
//
//
// Author(s) : Sebastien Loriot, Sylvain Pion
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/CGAL_Ipelet_base.h>
#include <CGAL/Delaunay_triangulation_2.h>
#include <CGAL/Triangulation_vertex_base_with_info_2.h>
#include <CGAL/Regular_triangulation_face_base_2.h>
#include <CGAL/Regular_triangulation_vertex_base_2.h>
#include <CGAL/Regular_triangulation_2.h>
#include "include/CGAL_ipelets/k_delaunay.h"
namespace CGAL_multi_delaunay{
typedef CGAL::Exact_predicates_inexact_constructions_kernel Kernel;
typedef CGAL::Delaunay_triangulation_2<Kernel> Delaunay;
typedef CGAL::Regular_triangulation_vertex_base_2<Kernel> Vb;
typedef CGAL::Triangulation_vertex_base_with_info_2<std::vector<Kernel::Point_2>,Kernel,Vb> VbI;
//~ typedef CGAL::Triangulation_vertex_base_with_info_2<std::list<Kernel::Point_2>,Kernel,Vb> VbI;
typedef CGAL::Regular_triangulation_face_base_2<Kernel > Fb;
typedef CGAL::Triangulation_data_structure_2<VbI,Fb> Tds;
typedef CGAL::Regular_triangulation_2<Kernel,Tds> RegularI;
typedef CGAL::Regular_triangulation_2<Kernel> Regular;
typedef Delaunay::Finite_vertices_iterator itVert;
typedef RegularI::Finite_vertices_iterator RitVert;
typedef Delaunay::Vertex_handle Vertex;
typedef RegularI::Vertex_handle VertexI;
typedef Delaunay::Vertex_circulator VCirculator;
typedef RegularI::Vertex_circulator RVCirculator;
const std::string sublabel[] = {
"Delaunay", "Delaunay 2", "Delaunay 3","Delaunay n-1", "Delaunay k", "Voronoi", "Voronoi 2", "Voronoi 3", "Voronoi n-1", "Voronoi k", "Help"
};
const std::string hlpmsg[] = {
"Generate k-th Delaunay triangulation and k-th dual Voronoi diagram. Note : k must be smaller than the number of input points."
};
class MdelaunayIpelet
: public CGAL::Ipelet_base<Kernel,11> {
public:
MdelaunayIpelet()
: CGAL::Ipelet_base<Kernel,11>("k order Delaunay",sublabel,hlpmsg){}
void protected_run(int);
};
void MdelaunayIpelet::protected_run(int fn)
{
Delaunay dt;
RegularI rti;
Regular rt;
//~ std::vector<Point_2> pt_list; I use instead pt_list
if (fn==10){
show_help(false);
return;
}
std::vector<Point_2> pt_list;
Iso_rectangle_2 bbox=read_active_objects( CGAL::dispatch_or_drop_output<Point_2>( std::back_inserter(pt_list) ) );
if (pt_list.empty()){
print_error_message("No mark selected");
return;
}
dt.insert(pt_list.begin(),pt_list.end());
switch(fn){
case 0://Classical Delauney
draw_in_ipe(dt);
break;
case 1:
CGAL_FALLTHROUGH;
case 6:
CGAL_FALLTHROUGH;
case 2:
CGAL_FALLTHROUGH;
case 7: //Delaunay and Voronoi for 2nd-3rd order
for (Delaunay::Finite_edges_iterator it=dt.finite_edges_begin();it!=dt.finite_edges_end();++it){
Point_2 pt0=it->first->vertex(Delaunay::cw(it->second))->point();
Point_2 pt1=it->first->vertex(Delaunay::ccw(it->second))->point();
VertexI vertI_cgal = rti.insert(Weighted_point_2(CGAL::midpoint(pt0,pt1),-CGAL::to_double(CGAL::squared_distance(pt0,pt1))/4.));
pt_list.clear();
pt_list.push_back(pt0);
pt_list.push_back(pt1);
vertI_cgal -> info() = pt_list;
}
if(fn==1){//Delauney 2 : just regular triangulation of all midpoints of delaunay segments with weight minus the squared length of the edge divided by 4
draw_in_ipe(rti);
break;
}
if(fn==2 || fn==7){ //Pour l'order 3
//CAN WE ITERATE OVER DELAUNEY TRIANGLES???
//WE MAY COUNT SEVERAL TIME SAME TRIANGLE WITH THE FOLLOWING METHOD
//iterate over adjacent point in the regular triangulation and compute a new wpoint for those having one common parent from delaunay
for (RegularI::Finite_edges_iterator it=rti.finite_edges_begin();it!=rti.finite_edges_end();++it){
Point_2 pt0_ori0=it->first->vertex(Delaunay::cw(it->second))->info().front();
Point_2 pt0_ori1=it->first->vertex(Delaunay::cw(it->second))->info().back();
Point_2 pt1_ori0=it->first->vertex(Delaunay::ccw(it->second))->info().front();
Point_2 pt1_ori1=it->first->vertex(Delaunay::ccw(it->second))->info().back();
if(CGAL::compare_xy(pt0_ori0,pt1_ori0)==CGAL::EQUAL || CGAL::compare_xy(pt0_ori1,pt1_ori0)==CGAL::EQUAL)
rt.insert(Weighted_point_2(CGAL::centroid(pt0_ori0,pt0_ori1,pt1_ori1),-CGAL::to_double(CGAL::squared_distance(pt0_ori0,pt0_ori1)+
CGAL::squared_distance(pt0_ori0,pt1_ori1)+CGAL::squared_distance(pt1_ori1,pt0_ori1))/9.));
else
if(CGAL::compare_xy(pt0_ori0,pt1_ori1)==CGAL::EQUAL || CGAL::compare_xy(pt0_ori1,pt1_ori1)==CGAL::EQUAL)
rt.insert(Weighted_point_2(CGAL::centroid(pt0_ori0,pt0_ori1,pt1_ori0),-CGAL::to_double(CGAL::squared_distance(pt0_ori0,pt0_ori1)+
CGAL::squared_distance(pt0_ori0,pt1_ori0)+CGAL::squared_distance(pt1_ori0,pt0_ori1))/9.));
}
if(fn==2){//Draw 3th Delauney
draw_in_ipe(rt);
break;
}
}
CGAL_FALLTHROUGH;
case 3://Delaunay and Voronoi of order n-1
CGAL_FALLTHROUGH;
case 8:
if(fn==3 ||fn==8){
int order = pt_list.size()-1;
double pt_x0 =0;//base to compute centroid of n-1 points
double pt_y0 =0;//base to compute centroid of n-1 points
double wt=0;//total weight : sum of all distances between two input points
for(std::vector<Point_2>::iterator it_pt = pt_list.begin();it_pt!=pt_list.end();++it_pt){
pt_x0 = pt_x0 + CGAL::to_double((*it_pt).x());
pt_y0 = pt_y0 + CGAL::to_double((*it_pt).y());
for(std::vector<Point_2>::iterator it_pt2 = it_pt+1;it_pt2!=pt_list.end();++it_pt2){
wt = wt + CGAL::to_double(CGAL::squared_distance(*it_pt,*it_pt2));
}
}
for(std::vector<Point_2>::iterator it_pt = pt_list.begin();it_pt!=pt_list.end();++it_pt){
double w = wt;
//compute centroid of the set of input points / *it_pt
double pt_x = pt_x0 - CGAL::to_double((*it_pt).x());
double pt_y = pt_y0 - CGAL::to_double((*it_pt).y());
//remove from w the sum of distance of input point from *it_pt
for(std::vector<Point_2>::iterator it_pt2 = pt_list.begin();it_pt2!=pt_list.end();++it_pt2){ //Weighted_point_2 equivalent
w = w - CGAL::to_double(CGAL::squared_distance(*it_pt,*it_pt2));
}
w = - w / (double) (order*order);
pt_x = pt_x / (double) order;
pt_y = pt_y / (double) order;
rt.insert(Weighted_point_2(Point_2(pt_x,pt_y),w));
}
if(fn==3){//Draw (n-1)th Delaunay
draw_in_ipe(rt);
break;
}
}
CGAL_FALLTHROUGH;
case 4:
CGAL_FALLTHROUGH;
case 9://k-th Delauney and Voronoi
if(fn==4 ||fn==9){
int order;
int ret_val;
std::tie(ret_val,order)=request_value_from_user<int>("Enter order");
if (ret_val < 0){
print_error_message("Incorrect value");
return;
}
int nb_pts = pt_list.size();
if(order<1 || order>=nb_pts){
print_error_message("Not a good order");
return;
}
k_delaunay<Kernel>(rt,pt_list,order);
if(fn==4){//Draw k-th delaunay
draw_in_ipe(rt);
break;
}
}
CGAL_FALLTHROUGH;
case 5://Draw Voronoi diagrams
if(fn==5) draw_dual_in_ipe(dt,bbox);
if(fn==6) draw_dual_in_ipe(rti,bbox);
if(fn==7 || fn==8 || fn==9) draw_dual_in_ipe(rt,bbox);
break;
}
}
}
CGAL_IPELET(CGAL_multi_delaunay::MdelaunayIpelet)
|