mirror of https://github.com/CGAL/cgal
remove weight calculator from the public API of fairing
it is still available internally
This commit is contained in:
parent
63a2f8af8e
commit
cde06c7088
|
|
@ -18,21 +18,33 @@ namespace Polygon_mesh_processing {
|
|||
|
||||
@tparam SparseLinearSolver a model of `SparseLinearAlgebraTraitsWithPreFactor_d`. If \ref thirdpartyEigen "Eigen" 3.2 (or greater) is available
|
||||
and `CGAL_EIGEN3_ENABLED` is defined, then an overload of `Eigen_solver_traits` is provided as default parameter.
|
||||
@tparam WeightCalculator a model of `FairWeightCalculator` and can be omitted to use default Cotangent weights
|
||||
@tparam PolygonMesh must be a model of `MutableFaceGraph`
|
||||
@tparam InputIterator iterator over input vertices
|
||||
|
||||
@param pmesh polygon mesh to be faired
|
||||
@param vertex_begin first iterator of the range of vertices
|
||||
@param vertex_end past-the-end iterator of the range of vertices
|
||||
@param weight_calculator function object to calculate weights, default to Cotangent weights and can be omitted
|
||||
@param continuity tangential continuity, default to `FAIRING_C_1` and can be omitted
|
||||
|
||||
@return `true` if fairing is successful, otherwise no vertex position is changed
|
||||
|
||||
@todo accuracy of solvers are not good, for example when there is no boundary condition pre_factor should fail, but it does not.
|
||||
\todo WeightCalculator should be a property map
|
||||
*/
|
||||
template<class SparseLinearSolver, class PolygonMesh, class InputIterator>
|
||||
bool fair(PolygonMesh& pmesh,
|
||||
InputIterator vb,
|
||||
InputIterator ve,
|
||||
Fairing_continuity continuity = FAIRING_C_1
|
||||
)
|
||||
{
|
||||
typedef CGAL::internal::Cotangent_weight_with_voronoi_area_fairing<PolygonMesh> Weight_calculator;
|
||||
return fair<SparseLinearSolver, Weight_calculator, PolygonMesh, InputIterator>
|
||||
(pmesh, vb, ve, Weight_calculator(pmesh), continuity);
|
||||
}
|
||||
|
||||
// use non-default weight calculator
|
||||
// WeightCalculator a model of `FairWeightCalculator`, can be omitted to use default Cotangent weights
|
||||
// weight_calculator a function object to calculate weights, defaults to Cotangent weights and can be omitted
|
||||
template<class SparseLinearSolver, class WeightCalculator, class PolygonMesh, class InputIterator>
|
||||
bool fair(PolygonMesh& pmesh,
|
||||
InputIterator vertex_begin,
|
||||
|
|
@ -60,19 +72,6 @@ namespace Polygon_mesh_processing {
|
|||
(pmesh, vb, ve, weight_calculator, continuity);
|
||||
}
|
||||
|
||||
//use default WeightCalculator
|
||||
template<class SparseLinearSolver, class PolygonMesh, class InputIterator>
|
||||
bool fair(PolygonMesh& pmesh,
|
||||
InputIterator vb,
|
||||
InputIterator ve,
|
||||
Fairing_continuity continuity = FAIRING_C_1
|
||||
)
|
||||
{
|
||||
typedef CGAL::internal::Cotangent_weight_with_voronoi_area_fairing<PolygonMesh> Weight_calculator;
|
||||
return fair<SparseLinearSolver, Weight_calculator, PolygonMesh, InputIterator>
|
||||
(pmesh, vb, ve, Weight_calculator(pmesh), continuity);
|
||||
}
|
||||
|
||||
//use default SparseLinearSolver and WeightCalculator
|
||||
template<class PolygonMesh, class InputIterator>
|
||||
bool fair(PolygonMesh& pmesh,
|
||||
|
|
|
|||
|
|
@ -116,7 +116,6 @@ namespace Polygon_mesh_processing {
|
|||
one can pass `CGAL::Default()` as `solver`.
|
||||
|
||||
@tparam SparseLinearSolver a model of `SparseLinearAlgebraTraitsWithPreFactor_d`
|
||||
@tparam WeightCalculator a model of `FairWeightCalculator`
|
||||
@tparam PolygonMesh must be model of `MutableFaceGraph`
|
||||
@tparam FaceOutputIterator iterator holding `boost::graph_traits<PolygonMesh>::%face_descriptor` for patch faces.
|
||||
@tparam VertexOutputIterator iterator holding `boost::graph_traits<PolygonMesh>::%vertex_descriptor` for patch vertices.
|
||||
|
|
@ -125,7 +124,6 @@ namespace Polygon_mesh_processing {
|
|||
@param border_halfedge a border halfedge incident to the hole
|
||||
@param face_out iterator over patch faces
|
||||
@param vertex_out iterator over patch vertices without including the boundary
|
||||
@param weight_calculator function object to calculate weights, default to Cotangent weights and can be omitted
|
||||
@param density_control_factor factor for density where larger values cause denser refinements
|
||||
@param continuity tangential continuity, default to `FAIRING_C_1` and can be omitted
|
||||
@param use_delaunay_triangulation if `true`, use the Delaunay triangulation face search space
|
||||
|
|
@ -138,8 +136,50 @@ namespace Polygon_mesh_processing {
|
|||
- @a vertex_out
|
||||
|
||||
\todo handle islands
|
||||
\todo WeightCalculator should be a property map
|
||||
*/
|
||||
template<typename SparseLinearSolver,
|
||||
typename PolygonMesh,
|
||||
typename FaceOutputIterator,
|
||||
typename VertexOutputIterator>
|
||||
boost::tuple<bool, FaceOutputIterator, VertexOutputIterator>
|
||||
triangulate_refine_and_fair_hole(PolygonMesh& pmesh,
|
||||
typename boost::graph_traits<PolygonMesh>::halfedge_descriptor border_halfedge,
|
||||
FaceOutputIterator face_out,
|
||||
VertexOutputIterator vertex_out,
|
||||
SparseLinearSolver
|
||||
#ifdef DOXYGEN_RUNNING
|
||||
solver,
|
||||
#else
|
||||
/* solver */,
|
||||
#endif
|
||||
double density_control_factor = std::sqrt(2.0),
|
||||
Fairing_continuity continuity = FAIRING_C_1,
|
||||
bool use_delaunay_triangulation = true)
|
||||
{
|
||||
bool use_dt3 =
|
||||
#ifdef CGAL_HOLE_FILLING_DO_NOT_USE_DT3
|
||||
false;
|
||||
#else
|
||||
use_delaunay_triangulation;
|
||||
#endif
|
||||
std::vector<typename boost::graph_traits<PolygonMesh>::vertex_descriptor> patch;
|
||||
|
||||
face_out = triangulate_and_refine_hole
|
||||
(pmesh, border_halfedge, face_out, std::back_inserter(patch),
|
||||
density_control_factor, use_dt3).first;
|
||||
|
||||
typedef internal::Fair_default_sparse_linear_solver::Solver Default_solver;
|
||||
typedef typename Default::Get<SparseLinearSolver, Default_solver>::type Solver;
|
||||
|
||||
bool fair_success = fair<Solver>(pmesh, patch.begin(), patch.end(), continuity);
|
||||
|
||||
vertex_out = std::copy(patch.begin(), patch.end(), vertex_out);
|
||||
return boost::make_tuple(fair_success, face_out, vertex_out);
|
||||
}
|
||||
|
||||
// use non-default weight calculator
|
||||
// WeightCalculator a model of `FairWeightCalculator`
|
||||
// weight_calculator function object to calculate weights, default to Cotangent weights and can be omitted
|
||||
template<class WeightCalculator,
|
||||
class SparseLinearSolver,
|
||||
class PolygonMesh,
|
||||
|
|
@ -201,10 +241,9 @@ namespace Polygon_mesh_processing {
|
|||
#else
|
||||
use_delaunay_triangulation;
|
||||
#endif
|
||||
CGAL::internal::Cotangent_weight_with_voronoi_area_fairing<PolygonMesh> wc(pmesh);
|
||||
|
||||
return triangulate_refine_and_fair_hole
|
||||
(pmesh, border_halfedge, face_out, vertex_out, wc, Default(),
|
||||
(pmesh, border_halfedge, face_out, vertex_out, Default(),
|
||||
density_control_factor, continuity, use_dt3);
|
||||
}
|
||||
|
||||
|
|
|
|||
|
|
@ -84,7 +84,7 @@ public slots:
|
|||
CGAL::Polygon_mesh_processing::fair(*selection_item->polyhedron(),
|
||||
selection_item->selected_vertices.begin(),
|
||||
selection_item->selected_vertices.end(),
|
||||
CGAL::internal::Cotangent_weight_with_voronoi_area_fairing<Polyhedron>(*selection_item->polyhedron()),
|
||||
// CGAL::internal::Cotangent_weight_with_voronoi_area_fairing<Polyhedron>(*selection_item->polyhedron()),
|
||||
continuity);
|
||||
selection_item->changed_with_poly_item();
|
||||
QApplication::restoreOverrideCursor();
|
||||
|
|
|
|||
Loading…
Reference in New Issue