Externalize the triangulation code for the facets in items

- Use an external class to triangulate the facets in Polyhedron_item, Surface_mesh_item, Polygon_soup_item and Polyhedron_selection_item.
This commit is contained in:
Maxime Gimeno 2016-06-03 15:54:50 +02:00
parent 128f042c7a
commit a3ed36ee25
6 changed files with 262 additions and 432 deletions

View File

@ -24,14 +24,8 @@
#include <CGAL/Polygon_mesh_processing/orientation.h>
#include <CGAL/Polygon_mesh_processing/repair.h>
#include <CGAL/Triangulation_vertex_base_with_info_2.h>
#include <CGAL/Triangulation_face_base_with_info_2.h>
#include <CGAL/Constrained_Delaunay_triangulation_2.h>
#include <CGAL/Constrained_triangulation_plus_2.h>
#include <CGAL/Triangulation_2_projection_traits_3.h>
#include <CGAL/Polygon_mesh_processing/compute_normal.h>
#include "triangulate_primitive.h"
struct Scene_polygon_soup_item_priv{
typedef Polygon_soup::Polygons::const_iterator Polygons_iterator;
@ -242,23 +236,7 @@ Scene_polygon_soup_item_priv::initializeBuffers(CGAL::Three::Viewer_interface* v
typedef Polyhedron::Traits Traits;
typedef Polygon_soup::Polygon_3 Facet;
typedef CGAL::Triangulation_2_projection_traits_3<Traits> P_traits;
typedef Polyhedron::Halfedge_handle Halfedge_handle;
struct Face_info {
Polyhedron::Halfedge_handle e[3];
bool is_external;
};
typedef CGAL::Triangulation_vertex_base_with_info_2<Halfedge_handle,
P_traits> Vb;
typedef CGAL::Triangulation_face_base_with_info_2<Face_info,
P_traits> Fb1;
typedef CGAL::Constrained_triangulation_face_base_2<P_traits, Fb1> Fb;
typedef CGAL::Triangulation_data_structure_2<Vb,Fb> TDS;
typedef CGAL::Exact_predicates_tag Itag;
typedef CGAL::Constrained_Delaunay_triangulation_2<P_traits,
TDS,
Itag> CDTbase;
typedef CGAL::Constrained_triangulation_plus_2<CDTbase> CDT;
void
Scene_polygon_soup_item_priv::triangulate_polygon(Polygons_iterator pit, int polygon_id) const
{
@ -268,8 +246,26 @@ Scene_polygon_soup_item_priv::triangulate_polygon(Polygons_iterator pit, int pol
const Point_3& pc = soup->points[pit->at(2)];
Traits::Vector_3 normal = CGAL::cross_product(pb-pa, pc -pa);
normal = normal / std::sqrt(normal * normal);
typedef FacetTriangulator<Polyhedron, Kernel, std::size_t> FT;
P_traits cdt_traits(normal);
double diagonal;
if(item->diagonalBbox() != std::numeric_limits<double>::infinity())
diagonal = item->diagonalBbox();
else
diagonal = 0.0;
std::size_t it = 0;
std::size_t it_end =pit->size();
std::vector<FT::PointAndId> pointIds;
do {
FT::PointAndId pointId;
pointId.point = soup->points[pit->at(it)];
pointId.id = pit->at(it);
pointIds.push_back(pointId);
} while( ++it != it_end );
FT triangulation(pointIds,normal,diagonal);
/* P_traits cdt_traits(normal);
CDT cdt(cdt_traits);
//A map used to associate the vertices in the triangulation to the points
@ -320,16 +316,14 @@ Scene_polygon_soup_item_priv::triangulate_polygon(Polygons_iterator pit, int pol
}
}
}
*/
//iterates on the internal faces to add the vertices to the positions
//and the normals to the appropriate vectors
int count =0;
for(CDT::Finite_faces_iterator
ffit = cdt.finite_faces_begin(),
end = cdt.finite_faces_end();
for(typename FT::CDT::Finite_faces_iterator
ffit = triangulation.cdt->finite_faces_begin(),
end = triangulation.cdt->finite_faces_end();
ffit != end; ++ffit)
{
count ++;
if(ffit->info().is_external)
continue;
@ -372,7 +366,7 @@ Scene_polygon_soup_item_priv::triangulate_polygon(Polygons_iterator pit, int pol
}
if(!soup->vcolors.empty())
{
CGAL::Color vcolor = soup->vcolors[p2p[ffit->vertex(i)]];
CGAL::Color vcolor = soup->vcolors[triangulation.v2v[ffit->vertex(i)]];
v_colors.push_back((float)vcolor.red()/255);
v_colors.push_back((float)vcolor.green()/255);
v_colors.push_back((float)vcolor.blue()/255);

View File

@ -8,13 +8,8 @@
#include <CGAL/AABB_tree.h>
#include <CGAL/AABB_traits.h>
#include <CGAL/Triangulation_vertex_base_with_info_2.h>
#include <CGAL/Triangulation_face_base_with_info_2.h>
#include <CGAL/Polyhedron_items_with_id_3.h>
#include <CGAL/Constrained_Delaunay_triangulation_2.h>
#include <CGAL/Triangulation_2_projection_traits_3.h>
#include <CGAL/Polygon_mesh_processing/compute_normal.h>
#include <CGAL/boost/graph/graph_traits_Polyhedron_3.h>
#include <CGAL/Polygon_mesh_processing/connected_components.h>
#include <CGAL/Polygon_mesh_processing/measure.h>
#include <CGAL/Polygon_mesh_processing/self_intersections.h>
@ -35,29 +30,11 @@
#include <boost/foreach.hpp>
#include <boost/container/flat_map.hpp>
#include "triangulate_primitive.h"
namespace PMP = CGAL::Polygon_mesh_processing;
typedef Polyhedron::Traits Traits;
typedef Polyhedron::Facet Facet;
typedef CGAL::Triangulation_2_projection_traits_3<Traits> P_traits;
typedef Polyhedron::Halfedge_handle Halfedge_handle;
struct Face_info {
Polyhedron::Halfedge_handle e[3];
bool is_external;
};
typedef CGAL::Triangulation_vertex_base_with_info_2<Halfedge_handle,
P_traits> Vb;
typedef CGAL::Triangulation_face_base_with_info_2<Face_info,
P_traits> Fb1;
typedef CGAL::Constrained_triangulation_face_base_2<P_traits, Fb1> Fb;
typedef CGAL::Triangulation_data_structure_2<Vb,Fb> TDS;
typedef CGAL::Exact_predicates_tag Itag;
typedef CGAL::Constrained_Delaunay_triangulation_2<P_traits,
TDS,
Itag> CDT;
//Make sure all the facets are triangles
typedef Polyhedron::Traits Kernel;
typedef Kernel::Point_3 Point;
typedef Kernel::Vector_3 Vector;
@ -153,6 +130,8 @@ struct Scene_polyhedron_item_priv{
delete targeted_id;
}
void* get_aabb_tree();
QList<Kernel::Triangle_3> triangulate_primitive(Polyhedron::Facet_iterator fit,
Traits::Vector_3 normal);
Color_vector colors_;
bool show_only_feature_edges_m;
bool show_feature_edges_m;
@ -210,9 +189,10 @@ struct Scene_polyhedron_item_priv{
const char* aabb_property_name = "Scene_polyhedron_item aabb tree";
QList<Kernel::Triangle_3> triangulate_primitive(Polyhedron::Facet_iterator fit,
QList<Kernel::Triangle_3> Scene_polyhedron_item_priv::triangulate_primitive(Polyhedron::Facet_iterator fit,
Traits::Vector_3 normal)
{
typedef FacetTriangulator<Polyhedron, Polyhedron::Traits, typename boost::graph_traits<Polyhedron>::vertex_descriptor> FT;
//The output list
QList<Kernel::Triangle_3> res;
//check if normal contains NaN values
@ -221,58 +201,17 @@ QList<Kernel::Triangle_3> triangulate_primitive(Polyhedron::Facet_iterator fit,
qDebug()<<"Warning in triangulation of the selection item: normal contains NaN values and is not valid.";
return QList<Kernel::Triangle_3>();
}
P_traits cdt_traits(normal);
CDT cdt(cdt_traits);
Facet::Halfedge_around_facet_circulator
he_circ = fit->facet_begin(),
he_circ_end(he_circ);
// Iterates on the vector of facet handles
typedef boost::graph_traits<Polyhedron>::vertex_descriptor vertex_descriptor;
boost::container::flat_map<CDT::Vertex_handle, vertex_descriptor> v2v;
CDT::Vertex_handle previous, first;
do {
CDT::Vertex_handle vh = cdt.insert(he_circ->vertex()->point());
v2v.insert(std::make_pair(vh, he_circ->vertex()));
if(first == 0) {
first = vh;
}
vh->info() = he_circ;
if(previous != 0 && previous != vh) {
cdt.insert_constraint(previous, vh);
}
previous = vh;
} while( ++he_circ != he_circ_end );
cdt.insert_constraint(previous, first);
// sets mark is_external
for(CDT::All_faces_iterator
fit2 = cdt.all_faces_begin(),
end = cdt.all_faces_end();
fit2 != end; ++fit2)
{
fit2->info().is_external = false;
}
//check if the facet is external or internal
std::queue<CDT::Face_handle> face_queue;
face_queue.push(cdt.infinite_vertex()->face());
while(! face_queue.empty() ) {
CDT::Face_handle fh = face_queue.front();
face_queue.pop();
if(fh->info().is_external) continue;
fh->info().is_external = true;
for(int i = 0; i <3; ++i) {
if(!cdt.is_constrained(std::make_pair(fh, i)))
{
face_queue.push(fh->neighbor(i));
}
}
}
double diagonal;
if(item->diagonalBbox() != std::numeric_limits<double>::infinity())
diagonal = item->diagonalBbox();
else
diagonal = 0.0;
FT triangulation(fit,normal,poly,diagonal);
//iterates on the internal faces to add the vertices to the positions
//and the normals to the appropriate vectors
for(CDT::Finite_faces_iterator
ffit = cdt.finite_faces_begin(),
end = cdt.finite_faces_end();
for(typename FT::CDT::Finite_faces_iterator
ffit = triangulation.cdt->finite_faces_begin(),
end = triangulation.cdt->finite_faces_end();
ffit != end; ++ffit)
{
if(ffit->info().is_external)
@ -280,8 +219,8 @@ QList<Kernel::Triangle_3> triangulate_primitive(Polyhedron::Facet_iterator fit,
res << Kernel::Triangle_3(ffit->vertex(0)->point(),
ffit->vertex(1)->point(),
ffit->vertex(2)->point());
ffit->vertex(1)->point(),
ffit->vertex(2)->point());
}
return res;
@ -358,6 +297,16 @@ void push_back_xyz(const TypeWithXYZ& t,
vector.push_back(t.z());
}
template<typename TypeWithRGB, typename ContainerWithPushBack>
void push_back_rgb(const TypeWithRGB& t,
ContainerWithPushBack& vector)
{
vector.push_back(t.redF());
vector.push_back(t.greenF());
vector.push_back(t.blueF());
}
bool Scene_polyhedron_item_priv::isFacetConvex(Facet_iterator f, const Polyhedron::Traits::Vector_3& N) const
{
typedef Polyhedron::Traits::Vector_3 Vector;
@ -484,155 +433,45 @@ void Scene_polyhedron_item_priv::triangulate_convex_facet(Facet_iterator f,
template<typename VertexNormalPmap>
void
Scene_polyhedron_item_priv::triangulate_facet(Scene_polyhedron_item::Facet_iterator fit,
const Traits::Vector_3& N,
const Traits::Vector_3& normal,
const VertexNormalPmap& vnmap,
const bool colors_only) const
{
Traits::Vector_3 normal = N;
Traits::Orientation orientation;
bool normal_is_ok;
Facet::Halfedge_around_facet_circulator
he = fit->facet_begin(),
he_end(he);
//Check if the normal is good, else find another one.
do{
normal_is_ok = true;
//Initializes the facet orientation
Polyhedron::Traits::Point_3 S,T;
S = source(he, *poly)->point();
T = target(he, *poly)->point();
Vector V1(S,T);
S = source(he->next(), *poly)->point();
T = target(he->next(), *poly)->point();
Vector V2(S,T);
//check if normal contains NaN values
if (normal.x() != normal.x() || normal.y() != normal.y() || normal.z() != normal.z()
|| normal == Traits::Vector_3(0,0,0))
normal_is_ok = false;
if(normal_is_ok)
{
orientation = Polyhedron::Traits::Orientation_3()(V1, V2, normal);
if( orientation == CGAL::COPLANAR )
normal_is_ok = false;
}
//Checks if the normal is good : if the normal is null
// or if it is coplanar to the facet, we need another one.
if(!normal_is_ok)
{
normal = CGAL::cross_product(V1,V2);
}
}while( ++he != he_end && !normal_is_ok);
//if no good normal can be found, stop here.
if (!normal_is_ok)
{
qDebug()<<"Warning : normal is not valid. Facet not displayed";
return;
}
P_traits cdt_traits(normal);
CDT cdt(cdt_traits);
Facet::Halfedge_around_facet_circulator
he_circ = fit->facet_begin(),
he_circ_end(he_circ);
// Iterates on the vector of facet handles
typedef boost::graph_traits<Polyhedron>::vertex_descriptor vertex_descriptor;
boost::container::flat_map<CDT::Vertex_handle, vertex_descriptor> v2v;
CDT::Vertex_handle previous, first, last_inserted;
/*Compute a reasonable precision level used to decide
* if two consecutive points in a facet can be estimated
* equal.*/
double min_sq_dist = CGAL::square(0.0001*item->DiagonalBbox());
// Iterate the points of the facet and decide if they must be inserted in the CDT
do {
/* If this is the first vertex, we insert it. Else, it must be reasonnably
* far from the previous vertex to be considered as a separated point.
* This is done to avoid problems with exact/inexact behavior, that causes
* the CDT's intersection calculation to crash.
*/
double sq_dist =
(first == 0) ? DBL_MAX :
(std::min)(
CGAL::squared_distance(previous->point(), he_circ->vertex()->point()),
CGAL::squared_distance(first->point(), he_circ->vertex()->point())
);
if( sq_dist < min_sq_dist)
continue;
CDT::Vertex_handle vh = cdt.insert(he_circ->vertex()->point());
v2v.insert(std::make_pair(vh, he_circ->vertex()));
if(first == 0) {
first = vh;
}
vh->info() = he_circ;
if(previous != 0 && previous != vh) {
cdt.insert_constraint(previous, vh);
last_inserted = previous;
}
previous = vh;
} while( ++he_circ != he_circ_end );
double sq_dist = CGAL::squared_distance(previous->point(), first->point());
if(sq_dist > min_sq_dist)
cdt.insert_constraint(previous, first);
typedef FacetTriangulator<Polyhedron, Polyhedron::Traits, typename boost::graph_traits<Polyhedron>::vertex_descriptor> FT;
double diagonal;
if(item->diagonalBbox() != std::numeric_limits<double>::infinity())
diagonal = item->diagonalBbox();
else
cdt.insert_constraint(last_inserted, first);
// sets mark is_external
for(CDT::All_faces_iterator
fit2 = cdt.all_faces_begin(),
end = cdt.all_faces_end();
fit2 != end; ++fit2)
{
fit2->info().is_external = false;
}
//check if the facet is external or internal
std::queue<CDT::Face_handle> face_queue;
face_queue.push(cdt.infinite_vertex()->face());
while(! face_queue.empty() ) {
CDT::Face_handle fh = face_queue.front();
face_queue.pop();
if(fh->info().is_external) continue;
fh->info().is_external = true;
for(int i = 0; i <3; ++i) {
if(!cdt.is_constrained(std::make_pair(fh, i)))
{
face_queue.push(fh->neighbor(i));
}
}
}
if(cdt.dimension() != 2 )
diagonal = 0.0;
FT triangulation(fit,normal,poly,diagonal);
if(triangulation.cdt->dimension() != 2 )
{
qDebug()<<"Warning : cdt not right. Facet not displayed";
return;
}
//iterates on the internal faces to add the vertices to the positions
//and the normals to the appropriate vectors
const int this_patch_id = fit->patch_id();
for(CDT::Finite_faces_iterator
ffit = cdt.finite_faces_begin(),
end = cdt.finite_faces_end();
for(typename FT::CDT::Finite_faces_iterator
ffit = triangulation.cdt->finite_faces_begin(),
end = triangulation.cdt->finite_faces_end();
ffit != end; ++ffit)
{
if(ffit->info().is_external)
continue;
if (item->isItemMulticolor())
{
for (int i = 0; i<3; ++i)
{
color_facets.push_back(colors_[this_patch_id-m_min_patch_id].redF());
color_facets.push_back(colors_[this_patch_id-m_min_patch_id].greenF());
color_facets.push_back(colors_[this_patch_id-m_min_patch_id].blueF());
color_facets.push_back(colors_[this_patch_id-m_min_patch_id].redF());
color_facets.push_back(colors_[this_patch_id-m_min_patch_id].greenF());
color_facets.push_back(colors_[this_patch_id-m_min_patch_id].blueF());
}
if (is_multicolor)
{
for (int i = 0; i<3; ++i)
{
push_back_rgb(colors_[this_patch_id-m_min_patch_id], color_facets);
}
}
if (colors_only)
continue;
continue;
push_back_xyz(ffit->vertex(0)->point(), positions_facets);
positions_facets.push_back(1.0);
@ -647,13 +486,13 @@ Scene_polyhedron_item_priv::triangulate_facet(Scene_polyhedron_item::Facet_itera
push_back_xyz(normal, normals_flat);
push_back_xyz(normal, normals_flat);
Traits::Vector_3 ng = get(vnmap, v2v[ffit->vertex(0)]);
typename Traits::Vector_3 ng = get(vnmap, triangulation.v2v[ffit->vertex(0)]);
push_back_xyz(ng, normals_gouraud);
ng = get(vnmap, v2v[ffit->vertex(1)]);
ng = get(vnmap, triangulation.v2v[ffit->vertex(1)]);
push_back_xyz(ng, normals_gouraud);
ng = get(vnmap, v2v[ffit->vertex(2)]);
ng = get(vnmap, triangulation.v2v[ffit->vertex(2)]);
push_back_xyz(ng, normals_gouraud);
}
}
@ -952,7 +791,7 @@ Scene_polyhedron_item_priv::compute_normals_and_vertices(const bool colors_only)
}
else
{
triangulate_facet(f, nf, nv_pmap, colors_only);
this->triangulate_facet(f, nf, nv_pmap, colors_only);
}
}

View File

@ -1,17 +1,10 @@
#include "Scene_polyhedron_selection_item.h"
#include <CGAL/Polygon_mesh_processing/compute_normal.h>
#include <CGAL/Triangulation_vertex_base_with_info_2.h>
#include <CGAL/Triangulation_face_base_with_info_2.h>
#include <CGAL/Constrained_Delaunay_triangulation_2.h>
#include <CGAL/Constrained_triangulation_plus_2.h>
#include <CGAL/Triangulation_2_projection_traits_3.h>
#include <CGAL/boost/graph/graph_traits_Polyhedron_3.h>
#include <boost/container/flat_map.hpp>
#include <CGAL/boost/graph/dijkstra_shortest_paths.h>
#include <CGAL/property_map.h>
#include <functional>
#include "triangulate_primitive.h"
struct Scene_polyhedron_selection_item_priv{
typedef Polyhedron::Vertex_handle Vertex_handle;
@ -39,8 +32,7 @@ struct Scene_polyhedron_selection_item_priv{
const Selection_set_vertex& p_sel_vertex, const Selection_set_facet &p_sel_facet, const Selection_set_edge &p_sel_edges) const;
void compute_temp_elements() const;
template<typename FaceNormalPmap>
void triangulate_facet(Facet_handle, const FaceNormalPmap&,
void triangulate_facet(Facet_handle, Kernel::Vector_3 normal,
std::vector<float> &p_facets,std::vector<float> &p_normals) const;
void tempInstructions(QString s1, QString s2);
@ -303,100 +295,30 @@ void push_back_xyz(const TypeWithXYZ& t,
typedef Polyhedron::Traits Traits;
typedef Polyhedron::Facet Facet;
typedef CGAL::Triangulation_2_projection_traits_3<Traits> P_traits;
typedef Polyhedron::Halfedge_handle Halfedge_handle;
struct Face_info {
Polyhedron::Halfedge_handle e[3];
bool is_external;
};
typedef CGAL::Triangulation_vertex_base_with_info_2<Halfedge_handle,
P_traits> Vb;
typedef CGAL::Triangulation_face_base_with_info_2<Face_info,
P_traits> Fb1;
typedef CGAL::Constrained_triangulation_face_base_2<P_traits, Fb1> Fb;
typedef CGAL::Triangulation_data_structure_2<Vb,Fb> TDS;
typedef CGAL::Exact_predicates_tag Itag;
typedef CGAL::Constrained_Delaunay_triangulation_2<P_traits,
TDS,
Itag> CDTbase;
typedef CGAL::Constrained_triangulation_plus_2<CDTbase> CDT;
//Make sure all the facets are triangles
typedef Polyhedron::Traits Kernel;
typedef Kernel::Point_3 Point;
typedef Kernel::Vector_3 Vector;
typedef Polyhedron::Facet_iterator Facet_iterator;
typedef Polyhedron::Halfedge_around_facet_circulator HF_circulator;
typedef Traits::Point_3 Point;
typedef Traits::Vector_3 Vector;
typedef Polyhedron::Halfedge_handle Halfedge_handle;
typedef boost::graph_traits<Polyhedron>::face_descriptor face_descriptor;
typedef boost::graph_traits<Polyhedron>::vertex_descriptor vertex_descriptor;
template<typename FaceNormalPmap>
void
Scene_polyhedron_selection_item_priv::triangulate_facet(Facet_handle fit,const FaceNormalPmap& fnmap,
Scene_polyhedron_selection_item_priv::triangulate_facet(Facet_handle fit,const Vector normal,
std::vector<float> &p_facets,std::vector<float> &p_normals ) const
{
//Computes the normal of the facet
Traits::Vector_3 normal = get(fnmap, fit);
//check if normal contains NaN values
if (normal.x() != normal.x() || normal.y() != normal.y() || normal.z() != normal.z())
{
qDebug()<<"Warning : normal is not valid. Facet not displayed";
return;
}
P_traits cdt_traits(normal);
CDT cdt(cdt_traits);
Facet::Halfedge_around_facet_circulator
he_circ = fit->facet_begin(),
he_circ_end(he_circ);
// Iterates on the vector of facet handles
typedef boost::graph_traits<Polyhedron>::vertex_descriptor vertex_descriptor;
boost::container::flat_map<CDT::Vertex_handle, vertex_descriptor> v2v;
CDT::Vertex_handle previous, first;
do {
CDT::Vertex_handle vh = cdt.insert(he_circ->vertex()->point());
v2v.insert(std::make_pair(vh, he_circ->vertex()));
if(first == 0) {
first = vh;
}
vh->info() = he_circ;
if(previous != 0 && previous != vh) {
cdt.insert_constraint(previous, vh);
}
previous = vh;
} while( ++he_circ != he_circ_end );
cdt.insert_constraint(previous, first);
// sets mark is_external
for(CDT::All_faces_iterator
fit2 = cdt.all_faces_begin(),
end = cdt.all_faces_end();
fit2 != end; ++fit2)
{
fit2->info().is_external = false;
}
//check if the facet is external or internal
std::queue<CDT::Face_handle> face_queue;
face_queue.push(cdt.infinite_vertex()->face());
while(! face_queue.empty() ) {
CDT::Face_handle fh = face_queue.front();
face_queue.pop();
if(fh->info().is_external) continue;
fh->info().is_external = true;
for(int i = 0; i <3; ++i) {
if(!cdt.is_constrained(std::make_pair(fh, i)))
{
face_queue.push(fh->neighbor(i));
}
}
}
typedef FacetTriangulator<Polyhedron, Polyhedron::Traits, typename boost::graph_traits<Polyhedron>::vertex_descriptor> FT;
double diagonal;
if(item->poly_item->diagonalBbox() != std::numeric_limits<double>::infinity())
diagonal = item->poly_item->diagonalBbox();
else
diagonal = 0.0;
FT triangulation(fit,normal,poly,diagonal);
//iterates on the internal faces to add the vertices to the positions
//and the normals to the appropriate vectors
for(CDT::Finite_faces_iterator
ffit = cdt.finite_faces_begin(),
end = cdt.finite_faces_end();
for(typename FT::CDT::Finite_faces_iterator
ffit = triangulation.cdt->finite_faces_begin(),
end = triangulation.cdt->finite_faces_end();
ffit != end; ++ffit)
{
if(ffit->info().is_external)
@ -438,32 +360,28 @@ if(!poly)
Facet_handle f = (*it);
if (f == boost::graph_traits<Polyhedron>::null_face())
continue;
Vector nf = get(nf_pmap, f);
if(is_triangle(f->halfedge(),*poly))
{
const Kernel::Vector_3 n =
CGAL::Polygon_mesh_processing::compute_face_normal(f, *item->poly_item->polyhedron());
p_normals.push_back(nf.x());
p_normals.push_back(nf.y());
p_normals.push_back(nf.z());
p_normals.push_back(n.x());
p_normals.push_back(n.y());
p_normals.push_back(n.z());
p_normals.push_back(nf.x());
p_normals.push_back(nf.y());
p_normals.push_back(nf.z());
p_normals.push_back(n.x());
p_normals.push_back(n.y());
p_normals.push_back(n.z());
p_normals.push_back(n.x());
p_normals.push_back(n.y());
p_normals.push_back(n.z());
p_normals.push_back(nf.x());
p_normals.push_back(nf.y());
p_normals.push_back(nf.z());
Polyhedron::Halfedge_around_facet_circulator
he = f->facet_begin(),
cend = he;
CGAL_For_all(he,cend)
{
const Kernel::Point_3& p = he->vertex()->point();
const Point& p = he->vertex()->point();
p_facets.push_back(p.x());
p_facets.push_back(p.y());
p_facets.push_back(p.z());
@ -471,8 +389,6 @@ if(!poly)
}
else if (is_quad(f->halfedge(), *poly))
{
Vector nf = get(nf_pmap, f);
//1st half-quad
Point p0 = f->halfedge()->vertex()->point();
Point p1 = f->halfedge()->next()->vertex()->point();
@ -501,7 +417,7 @@ if(!poly)
}
else
{
triangulate_facet(f, nf_pmap, p_facets, p_normals);
triangulate_facet(f, nf, p_facets, p_normals);
}
}
@ -509,8 +425,8 @@ if(!poly)
{
for(Selection_set_edge::iterator it = p_sel_edges.begin(); it != p_sel_edges.end(); ++it) {
const Kernel::Point_3& a = (it->halfedge())->vertex()->point();
const Kernel::Point_3& b = (it->halfedge())->opposite()->vertex()->point();
const Point& a = (it->halfedge())->vertex()->point();
const Point& b = (it->halfedge())->opposite()->vertex()->point();
p_lines.push_back(a.x());
p_lines.push_back(a.y());
p_lines.push_back(a.z());
@ -528,7 +444,7 @@ if(!poly)
end = p_sel_vertices.end();
it != end; ++it)
{
const Kernel::Point_3& p = (*it)->point();
const Point& p = (*it)->point();
p_points.push_back(p.x());
p_points.push_back(p.y());
p_points.push_back(p.z());
@ -554,7 +470,7 @@ void Scene_polyhedron_selection_item_priv::compute_temp_elements()const
end = item->fixed_vertices.end();
it != end; ++it)
{
const Kernel::Point_3& p = (*it)->point();
const Point& p = (*it)->point();
positions_fixed_points.push_back(p.x());
positions_fixed_points.push_back(p.y());
positions_fixed_points.push_back(p.z());

View File

@ -10,29 +10,12 @@
#include <CGAL/boost/graph/properties_Surface_mesh.h>
#include <CGAL/Surface_mesh/Surface_mesh.h>
#include <CGAL/Triangulation_vertex_base_with_info_2.h>
#include <CGAL/Triangulation_face_base_with_info_2.h>
#include <CGAL/Constrained_Delaunay_triangulation_2.h>
#include <CGAL/Constrained_triangulation_plus_2.h>
#include <CGAL/Triangulation_2_projection_traits_3.h>
#include <CGAL/Polygon_mesh_processing/compute_normal.h>
#include "triangulate_primitive.h"
typedef boost::graph_traits<Scene_surface_mesh_item::SMesh>::face_descriptor face_descriptor;
typedef boost::graph_traits<Scene_surface_mesh_item::SMesh>::halfedge_descriptor halfedge_descriptor;
typedef boost::graph_traits<Scene_surface_mesh_item::SMesh>::vertex_descriptor vertex_descriptor;
struct Face_info {
halfedge_descriptor e[3];
bool is_external;
};
typedef CGAL::Triangulation_2_projection_traits_3<Scene_surface_mesh_item::Kernel> P_traits;
typedef CGAL::Triangulation_vertex_base_with_info_2<halfedge_descriptor,P_traits> Vb;
typedef CGAL::Triangulation_face_base_with_info_2<Face_info, P_traits > Fb1;
typedef CGAL::Constrained_triangulation_face_base_2<P_traits, Fb1> Fb;
typedef CGAL::Triangulation_data_structure_2<Vb,Fb> TDS;
typedef CGAL::Exact_intersections_tag Itag;
typedef CGAL::Constrained_Delaunay_triangulation_2<P_traits, TDS, Itag> CDTbase;
typedef CGAL::Constrained_triangulation_plus_2<CDTbase> CDT;
struct Scene_surface_mesh_item_priv{
@ -560,58 +543,18 @@ Scene_surface_mesh_item_priv::triangulate_facet(face_descriptor fd,
qDebug()<<"Warning : normal is not valid. Facet not displayed";
return;
}
P_traits cdt_traits(normal);
CDT cdt(cdt_traits);
SMesh::Halfedge_around_face_circulator
he_circ(halfedge(fd,*smesh_), *smesh_),
he_circ_end(he_circ);
// Iterates on the vector of facet handles
boost::container::flat_map<CDT::Vertex_handle, vertex_descriptor> v2v;
CDT::Vertex_handle previous, first;
do {
CDT::Vertex_handle vh = cdt.insert(smesh_->point(source(*he_circ, *smesh_)));
if(index)
v2v.insert(std::make_pair(vh, source(*he_circ, *smesh_)));
if(first == 0) {
first = vh;
}
vh->info() = *he_circ;
if(previous != 0 && previous != vh) {
cdt.insert_constraint(previous, vh);
}
previous = vh;
} while( ++he_circ != he_circ_end );
cdt.insert_constraint(previous, first);
// sets mark is_external
for( CDT::All_faces_iterator
fit2 = cdt.all_faces_begin(),
end = cdt.all_faces_end();
fit2 != end; ++fit2)
{
fit2->info().is_external = false;
}
//check if the facet is external or internal
std::queue< CDT::Face_handle> face_queue;
face_queue.push(cdt.infinite_vertex()->face());
while(! face_queue.empty() ) {
CDT::Face_handle fh = face_queue.front();
face_queue.pop();
if(fh->info().is_external) continue;
fh->info().is_external = true;
for(int i = 0; i <3; ++i) {
if(!cdt.is_constrained(std::make_pair(fh, i)))
{
face_queue.push(fh->neighbor(i));
}
}
}
typedef FacetTriangulator<SMesh, Kernel, typename boost::graph_traits<SMesh>::vertex_descriptor> FT;
double diagonal;
if(item->diagonalBbox() != std::numeric_limits<double>::infinity())
diagonal = item->diagonalBbox();
else
diagonal = 0.0;
FT triangulation(fd,normal,smesh_,diagonal);
//iterates on the internal faces
for( CDT::Finite_faces_iterator
ffit = cdt.finite_faces_begin(),
end = cdt.finite_faces_end();
for( typename FT::CDT::Finite_faces_iterator
ffit = triangulation.cdt->finite_faces_begin(),
end = triangulation.cdt->finite_faces_end();
ffit != end; ++ffit)
{
if(ffit->info().is_external)
@ -640,9 +583,9 @@ Scene_surface_mesh_item_priv::triangulate_facet(face_descriptor fd,
//adds the indices to the appropriate vector
else
{
idx_data_.push_back((*im)[v2v[ffit->vertex(0)]]);
idx_data_.push_back((*im)[v2v[ffit->vertex(1)]]);
idx_data_.push_back((*im)[v2v[ffit->vertex(2)]]);
idx_data_.push_back((*im)[triangulation.v2v[ffit->vertex(0)]]);
idx_data_.push_back((*im)[triangulation.v2v[ffit->vertex(1)]]);
idx_data_.push_back((*im)[triangulation.v2v[ffit->vertex(2)]]);
}
}

View File

@ -0,0 +1,138 @@
#ifndef TRIANGULATE_PRIMITIVE
#define TRIANGULATE_PRIMITIVE
#include <CGAL/Triangulation_vertex_base_with_info_2.h>
#include <CGAL/Triangulation_face_base_with_info_2.h>
#include <CGAL/Constrained_Delaunay_triangulation_2.h>
#include <CGAL/Triangulation_2_projection_traits_3.h>
#include <CGAL/Three/Scene_item.h>
#include <CGAL/boost/graph/graph_traits_Polyhedron_3.h>
#include <CGAL/boost/graph/graph_traits_Surface_mesh.h>
#include <QColor>
#include <boost/container/flat_map.hpp>
//Make sure all the facets are triangles
template<class Mesh, typename Kernel, typename Index_type>
class FacetTriangulator
{
typedef Kernel Traits;
typedef typename boost::graph_traits<Mesh>::halfedge_descriptor halfedge_descriptor;
typedef typename Kernel::Vector_3 Vector;
typedef CGAL::Triangulation_2_projection_traits_3<Traits> P_traits;
typedef CGAL::Triangulation_vertex_base_with_info_2<halfedge_descriptor,
P_traits> Vb;
struct Face_info {
typename boost::graph_traits<Mesh>::halfedge_descriptor e[3];
bool is_external;
};
typedef CGAL::Triangulation_face_base_with_info_2<Face_info,
P_traits> Fb1;
typedef CGAL::Constrained_triangulation_face_base_2<P_traits, Fb1> Fb;
typedef CGAL::Triangulation_data_structure_2<Vb,Fb> TDS;
typedef CGAL::Exact_predicates_tag Itag;
public:
struct PointAndId {
typename Kernel::Point_3 point;
Index_type id;
};
typedef CGAL::Constrained_Delaunay_triangulation_2<P_traits,
TDS,
Itag> CDT;
CDT *cdt;
boost::container::flat_map<typename CDT::Vertex_handle, Index_type> v2v;
//Constructor
FacetTriangulator(typename boost::graph_traits<Mesh>::face_descriptor fd,
const Vector& normal,
Mesh *poly,
const double item_diag)
{
std::vector<PointAndId> idPoints;
BOOST_FOREACH(halfedge_descriptor he_circ, halfedges_around_face( halfedge(fd, *poly), *poly))
{
PointAndId idPoint;
idPoint.point = get(boost::vertex_point,*poly,source(he_circ, *poly));
idPoint.id = source(he_circ, *poly);
idPoints.push_back(idPoint);
}
triangulate(idPoints, normal, item_diag);
}
FacetTriangulator(std::vector<PointAndId > &idPoints,
const Vector& normal,
const double item_diag)
{
triangulate(idPoints, normal, item_diag);
}
private:
void triangulate( std::vector<PointAndId > &idPoints,
const Vector& normal,
const double item_diag )
{
P_traits cdt_traits(normal);
cdt = new CDT(cdt_traits);
typename CDT::Vertex_handle previous, first, last_inserted;
//Compute a reasonable precision level used to decide
//if two consecutive points in a facet can be estimated
//equal.
double min_sq_dist = CGAL::square(0.0001*item_diag);
// Iterate the points of the facet and decide if they must be inserted in the CDT
BOOST_FOREACH(PointAndId idPoint, idPoints)
{
typename CDT::Vertex_handle vh = cdt->insert(idPoint.point);
v2v[vh] = idPoint.id;
if(first == 0) {
first = vh;
}
if(previous != 0 && previous != vh) {
cdt->insert_constraint(previous, vh);
last_inserted = previous;
}
previous = vh;
}
double sq_dist = CGAL::squared_distance(previous->point(), first->point());
if(sq_dist > min_sq_dist)
cdt->insert_constraint(previous, first);
else
cdt->insert_constraint(last_inserted, first);
// sets mark is_external
for(typename CDT::All_faces_iterator
fit2 = cdt->all_faces_begin(),
end = cdt->all_faces_end();
fit2 != end; ++fit2)
{
fit2->info().is_external = false;
}
//check if the facet is external or internal
std::queue<typename CDT::Face_handle> face_queue;
face_queue.push(cdt->infinite_vertex()->face());
while(! face_queue.empty() ) {
typename CDT::Face_handle fh = face_queue.front();
face_queue.pop();
if(fh->info().is_external) continue;
fh->info().is_external = true;
for(int i = 0; i <3; ++i) {
if(!cdt->is_constrained(std::make_pair(fh, i)))
{
face_queue.push(fh->neighbor(i));
}
}
}
}
};
#endif // TRIANGULATE_PRIMITIVE

View File

@ -166,7 +166,7 @@ public:
is_bbox_computed = true;
return _bbox;
}
virtual double DiagonalBbox() const {
virtual double diagonalBbox() const {
if(!is_diag_bbox_computed)
compute_diag_bbox();
is_diag_bbox_computed = true;