mirror of https://github.com/CGAL/cgal
expanding from and evaluating on vertices instead of faces
This commit is contained in:
parent
6b985bfeb8
commit
12a627e23f
|
|
@ -17,7 +17,7 @@ typedef CGAL::Surface_mesh<EpicKernel::Point_3> Mesh;
|
|||
typedef boost::graph_traits<Mesh>::face_descriptor face_descriptor;
|
||||
typedef boost::graph_traits<Mesh>::vertex_descriptor vertex_descriptor;
|
||||
typedef std::unordered_map<face_descriptor, EpicKernel::FT> FaceMeasureMap_tag;
|
||||
typedef std::unordered_map<vertex_descriptor, EpicKernel::Vector_3> vertexVectorMap_tag;
|
||||
typedef std::unordered_map<vertex_descriptor, EpicKernel::FT> vertexVectorMap_tag;
|
||||
|
||||
|
||||
int main(int argc, char* argv[])
|
||||
|
|
@ -25,7 +25,7 @@ int main(int argc, char* argv[])
|
|||
Mesh g1;
|
||||
const std::string filename = (argc>1) ?
|
||||
argv[1] :
|
||||
CGAL::data_file_path("meshes/small_bunny.obj");
|
||||
CGAL::data_file_path("meshes/cylinder.off");
|
||||
|
||||
if(!CGAL::IO::read_polygon_mesh(filename, g1))
|
||||
{
|
||||
|
|
@ -33,8 +33,9 @@ int main(int argc, char* argv[])
|
|||
return EXIT_FAILURE;
|
||||
}
|
||||
|
||||
FaceMeasureMap_tag mean_curvature_init, gaussian_curvature_init;
|
||||
boost::associative_property_map<FaceMeasureMap_tag>
|
||||
|
||||
vertexVectorMap_tag mean_curvature_init, gaussian_curvature_init;
|
||||
boost::associative_property_map<vertexVectorMap_tag>
|
||||
mean_curvature_map(mean_curvature_init), gaussian_curvature_map(gaussian_curvature_init);
|
||||
|
||||
PMP::interpolated_corrected_mean_curvature(
|
||||
|
|
@ -46,7 +47,7 @@ int main(int argc, char* argv[])
|
|||
gaussian_curvature_map
|
||||
);
|
||||
|
||||
for (face_descriptor f: g1.faces())
|
||||
std::cout << f.idx() << ": HC = " << get(mean_curvature_map, f)
|
||||
<< ", GC = " << get(gaussian_curvature_map, f) << "\n";
|
||||
for (vertex_descriptor v : vertices(g1))
|
||||
std::cout << v.idx() << ": HC = " << get(mean_curvature_map, v)
|
||||
<< ", GC = " << get(gaussian_curvature_map, v) << "\n";
|
||||
}
|
||||
|
|
|
|||
|
|
@ -378,11 +378,12 @@ typename GT::FT face_in_ball_ratio_2(const std::vector<typename GT::Vector_3>& x
|
|||
return (r - d_min) / (d_max - d_min);
|
||||
}
|
||||
|
||||
template<typename PolygonMesh, typename FaceMeasureMap,
|
||||
template<typename PolygonMesh, typename FaceMeasureMap, typename VertexCurvatureMap,
|
||||
typename NamedParameters = parameters::Default_named_parameters>
|
||||
void expand_interpolated_corrected_measure_face(const PolygonMesh& pmesh,
|
||||
void expand_interpolated_corrected_measure_vertex(const PolygonMesh& pmesh,
|
||||
FaceMeasureMap fmm,
|
||||
const typename boost::graph_traits<PolygonMesh>::face_descriptor f,
|
||||
VertexCurvatureMap vcm,
|
||||
const typename boost::graph_traits<PolygonMesh>::vertex_descriptor v,
|
||||
const NamedParameters& np = parameters::default_values())
|
||||
{
|
||||
using parameters::choose_parameter;
|
||||
|
|
@ -406,30 +407,28 @@ template<typename PolygonMesh, typename FaceMeasureMap,
|
|||
std::queue<face_descriptor> bfs_q;
|
||||
std::unordered_set<face_descriptor> bfs_v;
|
||||
|
||||
//get face center c
|
||||
typename GT::Vector_3 c;
|
||||
for (vertex_descriptor v : vertices_around_face(halfedge(f, pmesh), pmesh))
|
||||
{
|
||||
typename GT::Point_3 p = get(vpm, v);
|
||||
c += typename GT::Vector_3(p.x(), p.y(), p.z());
|
||||
}
|
||||
c /= degree(f, pmesh);
|
||||
typename GT::Point_3 vp = get(vpm, v);
|
||||
typename GT::Vector_3 c = typename GT::Vector_3(vp.x(), vp.y(), vp.z());
|
||||
|
||||
typename GT::FT corrected_mui = 0;
|
||||
|
||||
bfs_q.push(f);
|
||||
bfs_v.insert(f);
|
||||
|
||||
for (face_descriptor f : faces_around_target(halfedge(v, pmesh), pmesh)) {
|
||||
if (f != boost::graph_traits<PolygonMesh>::null_face())
|
||||
{
|
||||
bfs_q.push(f);
|
||||
bfs_v.insert(f);
|
||||
}
|
||||
}
|
||||
while (!bfs_q.empty()) {
|
||||
face_descriptor fi = bfs_q.front();
|
||||
bfs_q.pop();
|
||||
|
||||
// looping over vertices in face to get point coordinates
|
||||
std::vector<typename GT::Vector_3> x;
|
||||
for (vertex_descriptor v : vertices_around_face(halfedge(fi, pmesh), pmesh))
|
||||
for (vertex_descriptor vi : vertices_around_face(halfedge(fi, pmesh), pmesh))
|
||||
{
|
||||
typename GT::Point_3 p = get(vpm, v);
|
||||
x.push_back(typename GT::Vector_3(p.x(), p.y(), p.z()));
|
||||
typename GT::Point_3 pi = get(vpm, vi);
|
||||
x.push_back(typename GT::Vector_3(pi.x(), pi.y(), pi.z()));
|
||||
}
|
||||
|
||||
const typename GT::FT f_ratio = face_in_ball_ratio_2<GT>(x, r, c);
|
||||
|
|
@ -448,7 +447,7 @@ template<typename PolygonMesh, typename FaceMeasureMap,
|
|||
}
|
||||
}
|
||||
|
||||
put(fmm, f, corrected_mui);
|
||||
put(vcm, v, corrected_mui);
|
||||
}
|
||||
|
||||
//template<typename PolygonMesh, typename FaceMeasureMap,
|
||||
|
|
@ -462,63 +461,75 @@ template<typename PolygonMesh, typename FaceMeasureMap,
|
|||
// expand_interpolated_corrected_measure_face(pmesh, fmm, f, np);
|
||||
//}
|
||||
|
||||
template<typename PolygonMesh, typename FaceCurvatureMap,
|
||||
template<typename PolygonMesh, typename VertexCurvatureMap,
|
||||
typename NamedParameters = parameters::Default_named_parameters>
|
||||
void interpolated_corrected_mean_curvature(const PolygonMesh& pmesh,
|
||||
FaceCurvatureMap fcm,
|
||||
VertexCurvatureMap vcm,
|
||||
const NamedParameters& np = parameters::default_values())
|
||||
{
|
||||
typedef typename GetGeomTraits<PolygonMesh, NamedParameters>::type GT;
|
||||
typedef typename boost::graph_traits<PolygonMesh>::face_descriptor face_descriptor;
|
||||
typedef std::unordered_map<face_descriptor, typename GT::FT> FaceMeasureMap_tag;
|
||||
typedef typename boost::graph_traits<PolygonMesh>::vertex_descriptor vertex_descriptor;
|
||||
typedef std::unordered_map<vertex_descriptor, typename GT::FT> VertexMeasureMap_tag;
|
||||
|
||||
FaceMeasureMap_tag mu0_init, mu1_init;
|
||||
boost::associative_property_map<FaceMeasureMap_tag>
|
||||
mu0_map(mu0_init), mu1_map(mu1_init);
|
||||
|
||||
VertexMeasureMap_tag mu0_expand_init, mu1_expand_init;
|
||||
boost::associative_property_map<VertexMeasureMap_tag>
|
||||
mu0_expand_map(mu0_expand_init), mu1_expand_map(mu1_expand_init);
|
||||
|
||||
interpolated_corrected_measure_mesh(pmesh, mu0_map, MU0_AREA_MEASURE);
|
||||
interpolated_corrected_measure_mesh(pmesh, mu1_map, MU1_MEAN_CURVATURE_MEASURE);
|
||||
|
||||
for (face_descriptor f : faces(pmesh))
|
||||
for (vertex_descriptor v : vertices(pmesh))
|
||||
{
|
||||
expand_interpolated_corrected_measure_face(pmesh, mu0_map, f, np);
|
||||
expand_interpolated_corrected_measure_face(pmesh, mu1_map, f, np);
|
||||
expand_interpolated_corrected_measure_vertex(pmesh, mu0_map, mu0_expand_map, v, np);
|
||||
expand_interpolated_corrected_measure_vertex(pmesh, mu1_map, mu1_expand_map, v, np);
|
||||
|
||||
const typename GT::FT f_mu0 = get(mu0_map, f);
|
||||
if (f_mu0 > 0.000001)
|
||||
put(fcm, f, 0.5 * get(mu1_map, f) / get(mu0_map, f));
|
||||
else
|
||||
put(fcm, f, 0);
|
||||
typename GT::FT v_mu0 = get(mu0_expand_map, v);
|
||||
if (v_mu0 > 0.000001)
|
||||
put(vcm, v, 0.5 * get(mu1_expand_map, v) / v_mu0);
|
||||
else
|
||||
put(vcm, v, 0);
|
||||
}
|
||||
}
|
||||
|
||||
template<typename PolygonMesh, typename FaceCurvatureMap,
|
||||
template<typename PolygonMesh, typename VertexCurvatureMap,
|
||||
typename NamedParameters = parameters::Default_named_parameters>
|
||||
void interpolated_corrected_gaussian_curvature(const PolygonMesh& pmesh,
|
||||
FaceCurvatureMap fcm,
|
||||
VertexCurvatureMap vcm,
|
||||
const NamedParameters& np = parameters::default_values())
|
||||
{
|
||||
typedef typename GetGeomTraits<PolygonMesh, NamedParameters>::type GT;
|
||||
typedef typename boost::graph_traits<PolygonMesh>::face_descriptor face_descriptor;
|
||||
typedef std::unordered_map<face_descriptor, typename GT::FT> FaceMeasureMap_tag;
|
||||
typedef typename boost::graph_traits<PolygonMesh>::vertex_descriptor vertex_descriptor;
|
||||
typedef std::unordered_map<vertex_descriptor, typename GT::FT> VertexMeasureMap_tag;
|
||||
|
||||
FaceMeasureMap_tag mu0_init, mu2_init;
|
||||
boost::associative_property_map<FaceMeasureMap_tag>
|
||||
mu0_map(mu0_init), mu2_map(mu2_init);
|
||||
|
||||
VertexMeasureMap_tag mu0_expand_init, mu2_expand_init;
|
||||
boost::associative_property_map<VertexMeasureMap_tag>
|
||||
mu0_expand_map(mu0_expand_init), mu2_expand_map(mu2_expand_init);
|
||||
|
||||
interpolated_corrected_measure_mesh(pmesh, mu0_map, MU0_AREA_MEASURE);
|
||||
interpolated_corrected_measure_mesh(pmesh, mu2_map, MU2_GAUSSIAN_CURVATURE_MEASURE);
|
||||
|
||||
for (face_descriptor f : faces(pmesh))
|
||||
for (vertex_descriptor v : vertices(pmesh))
|
||||
{
|
||||
expand_interpolated_corrected_measure_face(pmesh, mu0_map, f, np);
|
||||
expand_interpolated_corrected_measure_face(pmesh, mu2_map, f, np);
|
||||
expand_interpolated_corrected_measure_vertex(pmesh, mu0_map, mu0_expand_map, v, np);
|
||||
expand_interpolated_corrected_measure_vertex(pmesh, mu2_map, mu2_expand_map, v, np);
|
||||
|
||||
const typename GT::FT f_mu0 = get(mu0_map, f);
|
||||
if(f_mu0 > 0.000001)
|
||||
put(fcm, f, get(mu2_map, f) / f_mu0);
|
||||
typename GT::FT v_mu0 = get(mu0_expand_map, v);
|
||||
if(v_mu0 > 0.000001)
|
||||
put(vcm, v, get(mu2_expand_map, v) / v_mu0);
|
||||
else
|
||||
put(fcm, f, 0);
|
||||
put(vcm, v, 0);
|
||||
}
|
||||
}
|
||||
|
||||
|
|
|
|||
|
|
@ -611,11 +611,11 @@ private Q_SLOTS:
|
|||
break;
|
||||
case 4: // Interpolated Corrected Mean Curvature
|
||||
displayInterpolatedCurvatureMeasure(sm_item, PMP::MU1_MEAN_CURVATURE_MEASURE);
|
||||
sm_item->setRenderingMode(Flat);
|
||||
sm_item->setRenderingMode(Gouraud);
|
||||
break;
|
||||
case 5: // Interpolated Corrected Gaussian Curvature
|
||||
displayInterpolatedCurvatureMeasure(sm_item, PMP::MU2_GAUSSIAN_CURVATURE_MEASURE);
|
||||
sm_item->setRenderingMode(Flat);
|
||||
sm_item->setRenderingMode(Gouraud);
|
||||
break;
|
||||
default:
|
||||
if(dock_widget->propertyBox->currentText().contains("v:"))
|
||||
|
|
@ -652,11 +652,11 @@ private Q_SLOTS:
|
|||
if(does_exist)
|
||||
sm_item->face_graph()->remove_property_map(pmap);
|
||||
std::tie(pmap, does_exist) =
|
||||
sm_item->face_graph()->property_map<face_descriptor, double>("f:interpolated_corrected_mean_curvature");
|
||||
sm_item->face_graph()->property_map<face_descriptor, double>("v:interpolated_corrected_mean_curvature");
|
||||
if (does_exist)
|
||||
sm_item->face_graph()->remove_property_map(pmap);
|
||||
std::tie(pmap, does_exist) =
|
||||
sm_item->face_graph()->property_map<face_descriptor, double>("f:interpolated_corrected_gaussian_curvature");
|
||||
sm_item->face_graph()->property_map<face_descriptor, double>("v:interpolated_corrected_gaussian_curvature");
|
||||
if (does_exist)
|
||||
sm_item->face_graph()->remove_property_map(pmap);
|
||||
});
|
||||
|
|
@ -731,14 +731,14 @@ private Q_SLOTS:
|
|||
{
|
||||
smesh.remove_property_map(angles);
|
||||
}
|
||||
SMesh::Property_map<face_descriptor, double> mean_curvature;
|
||||
std::tie(mean_curvature, found) = smesh.property_map<face_descriptor, double>("f:interpolated_corrected_mean_curvature");
|
||||
SMesh::Property_map<vertex_descriptor, double> mean_curvature;
|
||||
std::tie(mean_curvature, found) = smesh.property_map<vertex_descriptor, double>("v:interpolated_corrected_mean_curvature");
|
||||
if (found)
|
||||
{
|
||||
smesh.remove_property_map(mean_curvature);
|
||||
}
|
||||
SMesh::Property_map<face_descriptor, double> gaussian_curvature;
|
||||
std::tie(gaussian_curvature, found) = smesh.property_map<face_descriptor, double>("f:interpolated_corrected_gaussian_curvature");
|
||||
SMesh::Property_map<vertex_descriptor, double> gaussian_curvature;
|
||||
std::tie(gaussian_curvature, found) = smesh.property_map<vertex_descriptor, double>("v:interpolated_corrected_gaussian_curvature");
|
||||
if (found)
|
||||
{
|
||||
smesh.remove_property_map(gaussian_curvature);
|
||||
|
|
@ -786,13 +786,13 @@ private Q_SLOTS:
|
|||
void displayInterpolatedCurvatureMeasure(Scene_surface_mesh_item* item, PMP::Curvature_measure_index mu_index)
|
||||
{
|
||||
std::string tied_string = (mu_index == PMP::MU1_MEAN_CURVATURE_MEASURE)?
|
||||
"f:interpolated_corrected_mean_curvature": "f:interpolated_corrected_gaussian_curvature";
|
||||
"v:interpolated_corrected_mean_curvature": "v:interpolated_corrected_gaussian_curvature";
|
||||
SMesh& smesh = *item->face_graph();
|
||||
//compute once and store the value per face
|
||||
bool non_init;
|
||||
SMesh::Property_map<face_descriptor, double> mu_i_map;
|
||||
SMesh::Property_map<vertex_descriptor, double> mu_i_map;
|
||||
std::tie(mu_i_map, non_init) =
|
||||
smesh.add_property_map<face_descriptor, double>(tied_string, 0);
|
||||
smesh.add_property_map<vertex_descriptor, double>(tied_string, 0);
|
||||
if (non_init)
|
||||
{
|
||||
if (mu_index == PMP::MU1_MEAN_CURVATURE_MEASURE)
|
||||
|
|
@ -802,18 +802,18 @@ private Q_SLOTS:
|
|||
|
||||
double res_min = ARBITRARY_DBL_MAX,
|
||||
res_max = -ARBITRARY_DBL_MAX;
|
||||
SMesh::Face_index min_index, max_index;
|
||||
for (SMesh::Face_index f : faces(smesh))
|
||||
SMesh::Vertex_index min_index, max_index;
|
||||
for (SMesh::Vertex_index v : vertices(smesh))
|
||||
{
|
||||
if (mu_i_map[f] > res_max)
|
||||
if (mu_i_map[v] > res_max)
|
||||
{
|
||||
res_max = mu_i_map[f];
|
||||
max_index = f;
|
||||
res_max = mu_i_map[v];
|
||||
max_index = v;
|
||||
}
|
||||
if (mu_i_map[f] < res_min)
|
||||
if (mu_i_map[v] < res_min)
|
||||
{
|
||||
res_min = mu_i_map[f];
|
||||
min_index = f;
|
||||
res_min = mu_i_map[v];
|
||||
min_index = v;
|
||||
}
|
||||
}
|
||||
switch (mu_index)
|
||||
|
|
@ -831,7 +831,7 @@ private Q_SLOTS:
|
|||
connect(item, &Scene_surface_mesh_item::itemChanged,
|
||||
this, &DisplayPropertyPlugin::resetProperty);
|
||||
}
|
||||
treat_sm_property<face_descriptor>(tied_string, item->face_graph());
|
||||
treat_sm_property<vertex_descriptor>(tied_string, item->face_graph());
|
||||
}
|
||||
|
||||
|
||||
|
|
@ -1199,7 +1199,7 @@ private Q_SLOTS:
|
|||
case 4:
|
||||
{
|
||||
::zoomToId(*item->face_graph(),
|
||||
QString("f%1").arg(mean_curvature_min[item].second),
|
||||
QString("v%1").arg(mean_curvature_min[item].second),
|
||||
getActiveViewer(),
|
||||
dummy_fd,
|
||||
dummy_p);
|
||||
|
|
@ -1208,7 +1208,7 @@ private Q_SLOTS:
|
|||
case 5:
|
||||
{
|
||||
::zoomToId(*item->face_graph(),
|
||||
QString("f%1").arg(gaussian_curvature_min[item].second),
|
||||
QString("v%1").arg(gaussian_curvature_min[item].second),
|
||||
getActiveViewer(),
|
||||
dummy_fd,
|
||||
dummy_p);
|
||||
|
|
@ -1252,7 +1252,7 @@ private Q_SLOTS:
|
|||
case 4:
|
||||
{
|
||||
::zoomToId(*item->face_graph(),
|
||||
QString("f%1").arg(mean_curvature_max[item].second),
|
||||
QString("v%1").arg(mean_curvature_max[item].second),
|
||||
getActiveViewer(),
|
||||
dummy_fd,
|
||||
dummy_p);
|
||||
|
|
@ -1261,7 +1261,7 @@ private Q_SLOTS:
|
|||
case 5:
|
||||
{
|
||||
::zoomToId(*item->face_graph(),
|
||||
QString("f%1").arg(gaussian_curvature_max[item].second),
|
||||
QString("v%1").arg(gaussian_curvature_max[item].second),
|
||||
getActiveViewer(),
|
||||
dummy_fd,
|
||||
dummy_p);
|
||||
|
|
@ -1557,11 +1557,11 @@ private:
|
|||
std::unordered_map<Scene_surface_mesh_item*, std::pair<double, SMesh::Face_index> > angles_min;
|
||||
std::unordered_map<Scene_surface_mesh_item*, std::pair<double, SMesh::Face_index> > angles_max;
|
||||
|
||||
std::unordered_map<Scene_surface_mesh_item*, std::pair<double, SMesh::Face_index> > mean_curvature_min;
|
||||
std::unordered_map<Scene_surface_mesh_item*, std::pair<double, SMesh::Face_index> > mean_curvature_max;
|
||||
std::unordered_map<Scene_surface_mesh_item*, std::pair<double, SMesh::Vertex_index> > mean_curvature_min;
|
||||
std::unordered_map<Scene_surface_mesh_item*, std::pair<double, SMesh::Vertex_index> > mean_curvature_max;
|
||||
|
||||
std::unordered_map<Scene_surface_mesh_item*, std::pair<double, SMesh::Face_index> > gaussian_curvature_min;
|
||||
std::unordered_map<Scene_surface_mesh_item*, std::pair<double, SMesh::Face_index> > gaussian_curvature_max;
|
||||
std::unordered_map<Scene_surface_mesh_item*, std::pair<double, SMesh::Vertex_index> > gaussian_curvature_min;
|
||||
std::unordered_map<Scene_surface_mesh_item*, std::pair<double, SMesh::Vertex_index> > gaussian_curvature_max;
|
||||
|
||||
std::unordered_map<Scene_surface_mesh_item*, Vertex_source_map> is_source;
|
||||
|
||||
|
|
|
|||
Loading…
Reference in New Issue