File: fair.h

package info (click to toggle)
cgal 6.1.1-1
  • links: PTS, VCS
  • area: main
  • in suites: forky
  • size: 144,952 kB
  • sloc: cpp: 811,597; ansic: 208,576; sh: 493; python: 411; makefile: 286; javascript: 174
file content (188 lines) | stat: -rw-r--r-- 8,222 bytes parent folder | download | duplicates (2)
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
// Copyright (c) 2015 GeometryFactory (France).
// All rights reserved.
//
// This file is part of CGAL (www.cgal.org).
//
// $URL: https://github.com/CGAL/cgal/blob/v6.1.1/Polygon_mesh_processing/include/CGAL/Polygon_mesh_processing/fair.h $
// $Id: include/CGAL/Polygon_mesh_processing/fair.h 08b27d3db14 $
// SPDX-License-Identifier: GPL-3.0-or-later OR LicenseRef-Commercial
//
//
// Author(s)     : Ilker O. Yaz

#ifndef CGAL_POLYGON_MESH_PROCESSING_FAIR_H
#define CGAL_POLYGON_MESH_PROCESSING_FAIR_H

#include <CGAL/license/Polygon_mesh_processing/meshing_hole_filling.h>

#include <CGAL/disable_warnings.h>

#include <CGAL/Polygon_mesh_processing/internal/fair_impl.h>
#include <CGAL/Named_function_parameters.h>
#include <CGAL/boost/graph/named_params_helper.h>
#include <CGAL/Weights/cotangent_weights.h>

#if defined(CGAL_EIGEN3_ENABLED)
#include <CGAL/Eigen_solver_traits.h>  // for sparse linear system solver
#endif

#include <type_traits>

namespace CGAL {

namespace Polygon_mesh_processing {

namespace internal {
  // use non-default weight calculator and non-default solver
  // 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<typename SparseLinearSolver,
           typename WeightCalculator,
           typename TriangleMesh,
           typename VertexRange,
           typename VertexPointMap>
bool fair(TriangleMesh& tmesh,
          const VertexRange& vertices,
          SparseLinearSolver solver,
          WeightCalculator weight_calculator,
          unsigned int continuity,
          VertexPointMap vpmap)
{
  CGAL::Polygon_mesh_processing::internal::Fair_Polyhedron_3
     <TriangleMesh, SparseLinearSolver, WeightCalculator, VertexPointMap>
     fair_functor(tmesh, vpmap, weight_calculator);
  return fair_functor.fair(vertices, solver, continuity);
}

} //end namespace internal

  /*!
  \ingroup PMP_meshing_grp

  @brief fairs a region on a triangle mesh.

  The points of the selected vertices are
  relocated to yield an as-smooth-as-possible surface patch,
  based on solving a linear bi-Laplacian system with boundary constraints,
  described in \cgalCite{Botsch2008OnLinearVariational}.
  The optional parameter `fairing_continuity` gives the ability to control the tangential
  continuity C<sup> n</sup> of the output mesh.

  The region described by `vertices` might contain multiple disconnected components.
  Note that the mesh connectivity is not altered in any way,
  only vertex locations get updated.

  Fairing might fail if fixed vertices, which are used as boundary conditions,
  do not suffice to solve constructed linear system.

  Note that if the vertex range to which fairing is applied contains all the vertices of the triangle mesh,
  fairing does not fail, but the mesh gets shrunk to `CGAL::ORIGIN`.

  @tparam TriangleMesh a model of `FaceGraph` and `MutableFaceGraph`
  @tparam VertexRange a range of vertex descriptors of `TriangleMesh`, model of `Range`.
          Its iterator type is `InputIterator`.
  @tparam NamedParameters a sequence of \ref bgl_namedparameters "Named Parameters"

  @param tmesh the triangle mesh with patches to be faired
  @param vertices the vertices of the patches to be faired (the positions of only those vertices will be changed)
  @param np an optional sequence of \ref bgl_namedparameters "Named Parameters" among the ones listed below

  \cgalNamedParamsBegin
    \cgalParamNBegin{vertex_point_map}
      \cgalParamDescription{a property map associating points to the vertices of `tmesh`}
      \cgalParamType{a class model of `ReadablePropertyMap` with `boost::graph_traits<TriangleMesh>::%vertex_descriptor`
                     as key type and `%Point_3` as value type}
      \cgalParamDefault{`boost::get(CGAL::vertex_point, tmesh)`}
      \cgalParamExtra{If this parameter is omitted, an internal property map for `CGAL::vertex_point_t`
                      should be available for the vertices of `tmesh`.}
    \cgalParamNEnd

    \cgalParamNBegin{fairing_continuity}
      \cgalParamDescription{A value controlling the tangential continuity of the output surface patch.
                            The possible values are 0, 1 and 2, referring to the  C<sup>0</sup>, C<sup>1</sup>
                            and C<sup>2</sup> continuity.}
      \cgalParamType{unsigned int}
      \cgalParamDefault{`1`}
      \cgalParamExtra{The larger `fairing_continuity` gets, the more fixed vertices are required.}
    \cgalParamNEnd

    \cgalParamNBegin{sparse_linear_solver}
      \cgalParamDescription{an instance of the sparse linear solver used for fairing}
      \cgalParamType{a class model of `SparseLinearAlgebraWithFactorTraits_d`}
      \cgalParamDefault{If \ref thirdpartyEigen "Eigen" 3.2 (or greater) is available and
                        `CGAL_EIGEN3_ENABLED` is defined, then the following overload of `Eigen_solver_traits`
                        is provided as default value:\n
                        `CGAL::Eigen_solver_traits<Eigen::SparseLU<CGAL::Eigen_sparse_matrix<double>::%EigenType, Eigen::COLAMDOrdering<int> > >`}
    \cgalParamNEnd
  \cgalNamedParamsEnd

  @return `true` if fairing is successful, otherwise no vertices are relocated.

  @pre `is_triangle_mesh(tmesh)`

  @warning This function involves linear algebra, that is computed using non-exact, floating-point arithmetic.

  @todo accuracy of solvers are not good, for example when there is no boundary condition pre_factor should fail, but it does not.
  */
  template<typename TriangleMesh,
           typename VertexRange,
           typename NamedParameters = parameters::Default_named_parameters>
  bool fair(TriangleMesh& tmesh,
            const VertexRange& vertices,
            const NamedParameters& np = parameters::default_values())
  {
    using parameters::get_parameter;
    using parameters::choose_parameter;

    CGAL_precondition(is_triangle_mesh(tmesh));

#if defined(CGAL_EIGEN3_ENABLED)
  #if EIGEN_VERSION_AT_LEAST(3,2,0)
    typedef CGAL::Eigen_solver_traits<Eigen::SparseLU<
      CGAL::Eigen_sparse_matrix<double>::EigenType, Eigen::COLAMDOrdering<int> >  >
      Default_solver;
  #else
    typedef bool Default_solver;//compilation should crash
      //if no solver is provided and Eigen version < 3.2
  #endif
#else
    typedef bool Default_solver;//compilation should crash
      //if no solver is provided and Eigen version < 3.2
#endif

#if defined(CGAL_EIGEN3_ENABLED)
    static_assert(
      (!std::is_same<typename GetSolver<NamedParameters, Default_solver>::type, bool>::value) || EIGEN_VERSION_AT_LEAST(3, 2, 0),
      "The function `fair` requires Eigen3 version 3.2 or later.");
#else
    static_assert(
      (!std::is_same<typename GetSolver<NamedParameters, Default_solver>::type, bool>::value),
      "The function `fair` requires Eigen3 version 3.2 or later.");
#endif

    typedef typename GetVertexPointMap < TriangleMesh, NamedParameters>::type VPMap;
    VPMap vpmap = choose_parameter(get_parameter(np, internal_np::vertex_point),
                                   get_property_map(vertex_point, tmesh));

    typedef typename GetGeomTraits < TriangleMesh, NamedParameters>::type GT;
    GT gt = choose_parameter<GT>(get_parameter(np, internal_np::geom_traits));

    // Cotangent_weight_with_voronoi_area_fairing has been changed to the version:
    // Secure_cotangent_weight_with_voronoi_area to avoid imprecisions from
    // the issue #4706 - https://github.com/CGAL/cgal/issues/4706.
    typedef CGAL::Weights::Secure_cotangent_weight_with_voronoi_area<TriangleMesh, VPMap, GT> Default_weight_calculator;

    return internal::fair(tmesh, vertices,
      choose_parameter<Default_solver>(get_parameter(np, internal_np::sparse_linear_solver)),
      choose_parameter(get_parameter(np, internal_np::weight_calculator), Default_weight_calculator(tmesh, vpmap, gt)),
      choose_parameter(get_parameter(np, internal_np::fairing_continuity), 1),
      vpmap);
  }

} // namespace Polygon_mesh_processing

} // namespace CGAL

#include <CGAL/enable_warnings.h>

#endif //CGAL_POLYGON_MESH_PROCESSING_FAIR_H