Point Cloud Library (PCL)  1.14.0-dev
poisson.h
1 /*
2  * Software License Agreement (BSD License)
3  *
4  * Copyright (c) 2010, Willow Garage, Inc.
5  * All rights reserved.
6  *
7  * Redistribution and use in source and binary forms, with or without
8  * modification, are permitted provided that the following conditions
9  * are met:
10  *
11  * * Redistributions of source code must retain the above copyright
12  * notice, this list of conditions and the following disclaimer.
13  * * Redistributions in binary form must reproduce the above
14  * copyright notice, this list of conditions and the following
15  * disclaimer in the documentation and/or other materials provided
16  * with the distribution.
17  * * Neither the name of Willow Garage, Inc. nor the names of its
18  * contributors may be used to endorse or promote products derived
19  * from this software without specific prior written permission.
20  *
21  * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
22  * "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
23  * LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
24  * FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
25  * COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
26  * INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
27  * BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
28  * LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
29  * CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
30  * LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
31  * ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
32  * POSSIBILITY OF SUCH DAMAGE.
33  *
34  * $Id$
35  *
36  */
37 
38 #pragma once
39 
40 #include <pcl/memory.h>
41 #include <pcl/pcl_macros.h>
42 #include <pcl/surface/reconstruction.h>
43 
44 namespace pcl
45 {
46  namespace poisson
47  {
48  class CoredVectorMeshData;
49  template <class Real> struct Point3D;
50  }
51 
52  /** \brief The Poisson surface reconstruction algorithm.
53  * \note Code adapted from Misha Kazhdan: http://www.cs.jhu.edu/~misha/Code/PoissonRecon/
54  * \note Based on the paper:
55  * * Michael Kazhdan, Matthew Bolitho, Hugues Hoppe, "Poisson surface reconstruction",
56  * SGP '06 Proceedings of the fourth Eurographics symposium on Geometry processing
57  * \author Alexandru-Eugen Ichim
58  * \ingroup surface
59  */
60  template<typename PointNT>
61  class Poisson : public SurfaceReconstruction<PointNT>
62  {
63  public:
64  using Ptr = shared_ptr<Poisson<PointNT> >;
65  using ConstPtr = shared_ptr<const Poisson<PointNT> >;
66 
69 
71 
73  using KdTreePtr = typename KdTree::Ptr;
74 
75  /** \brief Constructor that sets all the parameters to working default values. */
76  Poisson ();
77 
78  /** \brief Destructor. */
79  ~Poisson () override;
80 
81  /** \brief Create the surface.
82  * \param[out] output the resultant polygonal mesh
83  */
84  void
85  performReconstruction (pcl::PolygonMesh &output) override;
86 
87  /** \brief Create the surface.
88  * \param[out] points the vertex positions of the resulting mesh
89  * \param[out] polygons the connectivity of the resulting mesh
90  */
91  void
93  std::vector<pcl::Vertices> &polygons) override;
94 
95  /** \brief Set the maximum depth of the tree that will be used for surface reconstruction.
96  * \note Running at depth d corresponds to solving on a voxel grid whose resolution is no larger than
97  * 2^d x 2^d x 2^d. Note that since the reconstructor adapts the octree to the sampling density, the specified
98  * reconstruction depth is only an upper bound.
99  * \param[in] depth the depth parameter
100  */
101  inline void
102  setDepth (int depth) { depth_ = depth; }
103 
104  /** \brief Get the depth parameter */
105  inline int
106  getDepth () { return depth_; }
107 
108  inline void
109  setMinDepth (int min_depth) { min_depth_ = min_depth; }
110 
111  inline int
112  getMinDepth () { return min_depth_; }
113 
114  inline void
115  setPointWeight (float point_weight) { point_weight_ = point_weight; }
116 
117  inline float
118  getPointWeight () { return point_weight_; }
119 
120  /** \brief Set the ratio between the diameter of the cube used for reconstruction and the diameter of the
121  * samples' bounding cube.
122  * \param[in] scale the given parameter value
123  */
124  inline void
125  setScale (float scale) { scale_ = scale; }
126 
127  /** Get the ratio between the diameter of the cube used for reconstruction and the diameter of the
128  * samples' bounding cube.
129  */
130  inline float
131  getScale () { return scale_; }
132 
133  /** \brief Set the the depth at which a block Gauss-Seidel solver is used to solve the Laplacian equation
134  * \note Using this parameter helps reduce the memory overhead at the cost of a small increase in
135  * reconstruction time. (In practice, we have found that for reconstructions of depth 9 or higher a subdivide
136  * depth of 7 or 8 can greatly reduce the memory usage.)
137  * \param[in] solver_divide the given parameter value
138  */
139  inline void
140  setSolverDivide (int solver_divide) { solver_divide_ = solver_divide; }
141 
142  /** \brief Get the depth at which a block Gauss-Seidel solver is used to solve the Laplacian equation */
143  inline int
144  getSolverDivide () { return solver_divide_; }
145 
146  /** \brief Set the depth at which a block iso-surface extractor should be used to extract the iso-surface
147  * \note Using this parameter helps reduce the memory overhead at the cost of a small increase in extraction
148  * time. (In practice, we have found that for reconstructions of depth 9 or higher a subdivide depth of 7 or 8
149  * can greatly reduce the memory usage.)
150  * \param[in] iso_divide the given parameter value
151  */
152  inline void
153  setIsoDivide (int iso_divide) { iso_divide_ = iso_divide; }
154 
155  /** \brief Get the depth at which a block iso-surface extractor should be used to extract the iso-surface */
156  inline int
157  getIsoDivide () { return iso_divide_; }
158 
159  /** \brief Set the minimum number of sample points that should fall within an octree node as the octree
160  * construction is adapted to sampling density
161  * \note For noise-free samples, small values in the range [1.0 - 5.0] can be used. For more noisy samples,
162  * larger values in the range [15.0 - 20.0] may be needed to provide a smoother, noise-reduced, reconstruction.
163  * \param[in] samples_per_node the given parameter value
164  */
165  inline void
166  setSamplesPerNode (float samples_per_node) { samples_per_node_ = samples_per_node; }
167 
168  /** \brief Get the minimum number of sample points that should fall within an octree node as the octree
169  * construction is adapted to sampling density
170  */
171  inline float
172  getSamplesPerNode () { return samples_per_node_; }
173 
174  /** \brief Set the confidence flag
175  * \note Enabling this flag tells the reconstructor to use the size of the normals as confidence information.
176  * When the flag is not enabled, all normals are normalized to have unit-length prior to reconstruction.
177  * \param[in] confidence the given flag
178  */
179  inline void
180  setConfidence (bool confidence) { confidence_ = confidence; }
181 
182  /** \brief Get the confidence flag */
183  inline bool
184  getConfidence () { return confidence_; }
185 
186  /** \brief Enabling this flag tells the reconstructor to output a polygon mesh (rather than triangulating the
187  * results of Marching Cubes).
188  * \param[in] output_polygons the given flag
189  */
190  inline void
191  setOutputPolygons (bool output_polygons) { output_polygons_ = output_polygons; }
192 
193  /** \brief Get whether the algorithm outputs a polygon mesh or a triangle mesh */
194  inline bool
195  getOutputPolygons () { return output_polygons_; }
196 
197  /** \brief Set the degree parameter
198  * \param[in] degree the given degree
199  */
200  inline void
201  setDegree (int degree) { degree_ = degree; }
202 
203  /** \brief Get the degree parameter */
204  inline int
205  getDegree () { return degree_; }
206 
207  /** \brief Set the manifold flag.
208  * \note Enabling this flag tells the reconstructor to add the polygon barycenter when triangulating polygons
209  * with more than three vertices.
210  * \param[in] manifold the given flag
211  */
212  inline void
213  setManifold (bool manifold) { manifold_ = manifold; }
214 
215  /** \brief Get the manifold flag */
216  inline bool
217  getManifold () { return manifold_; }
218 
219  /** \brief Set the number of threads to use.
220  * \param[in] threads the number of threads
221  */
222  void
223  setThreads(int threads);
224 
225 
226  /** \brief Get the number of threads*/
227  inline int
229  {
230  return threads_;
231  }
232 
233  protected:
234  /** \brief Class get name method. */
235  std::string
236  getClassName () const override { return ("Poisson"); }
237 
238  private:
239  int depth_{8};
240  int min_depth_{5};
241  float point_weight_{4};
242  float scale_{1.1f};
243  int solver_divide_{8};
244  int iso_divide_{8};
245  float samples_per_node_{1.0};
246  bool confidence_{false};
247  bool output_polygons_{false};
248 
249  bool no_reset_samples_{false};
250  bool no_clip_tree_{false};
251  bool manifold_{true};
252 
253  int refine_{3};
254  int kernel_depth_{8};
255  int degree_{2};
256  bool non_adaptive_weights_{false};
257  bool show_residual_{false};
258  int min_iterations_{8};
259  float solver_accuracy_{1e-3f};
260  int threads_{1};
261 
262  template<int Degree> void
263  execute (poisson::CoredVectorMeshData &mesh,
264  poisson::Point3D<float> &translate,
265  float &scale);
266 
267  public:
269  };
270 }
271 
272 #ifdef PCL_NO_PRECOMPILE
273 #include <pcl/surface/impl/poisson.hpp>
274 #endif
KdTree represents the base spatial locator class for kd-tree implementations.
Definition: kdtree.h:56
shared_ptr< KdTree< PointT > > Ptr
Definition: kdtree.h:69
shared_ptr< PointCloud< PointT > > Ptr
Definition: point_cloud.h:413
The Poisson surface reconstruction algorithm.
Definition: poisson.h:62
void setThreads(int threads)
Set the number of threads to use.
Definition: poisson.hpp:71
std::string getClassName() const override
Class get name method.
Definition: poisson.h:236
int getThreads()
Get the number of threads.
Definition: poisson.h:228
void setSamplesPerNode(float samples_per_node)
Set the minimum number of sample points that should fall within an octree node as the octree construc...
Definition: poisson.h:166
float getPointWeight()
Definition: poisson.h:118
void setDegree(int degree)
Set the degree parameter.
Definition: poisson.h:201
Poisson()
Constructor that sets all the parameters to working default values.
typename pcl::PointCloud< PointNT >::Ptr PointCloudPtr
Definition: poisson.h:70
int getDepth()
Get the depth parameter.
Definition: poisson.h:106
void performReconstruction(pcl::PolygonMesh &output) override
Create the surface.
Definition: poisson.hpp:142
void setIsoDivide(int iso_divide)
Set the depth at which a block iso-surface extractor should be used to extract the iso-surface.
Definition: poisson.h:153
float getSamplesPerNode()
Get the minimum number of sample points that should fall within an octree node as the octree construc...
Definition: poisson.h:172
void setConfidence(bool confidence)
Set the confidence flag.
Definition: poisson.h:180
bool getConfidence()
Get the confidence flag.
Definition: poisson.h:184
void setPointWeight(float point_weight)
Definition: poisson.h:115
int getSolverDivide()
Get the depth at which a block Gauss-Seidel solver is used to solve the Laplacian equation.
Definition: poisson.h:144
int getDegree()
Get the degree parameter.
Definition: poisson.h:205
void setOutputPolygons(bool output_polygons)
Enabling this flag tells the reconstructor to output a polygon mesh (rather than triangulating the re...
Definition: poisson.h:191
void setDepth(int depth)
Set the maximum depth of the tree that will be used for surface reconstruction.
Definition: poisson.h:102
bool getOutputPolygons()
Get whether the algorithm outputs a polygon mesh or a triangle mesh.
Definition: poisson.h:195
float getScale()
Get the ratio between the diameter of the cube used for reconstruction and the diameter of the sample...
Definition: poisson.h:131
bool getManifold()
Get the manifold flag.
Definition: poisson.h:217
void setScale(float scale)
Set the ratio between the diameter of the cube used for reconstruction and the diameter of the sample...
Definition: poisson.h:125
typename KdTree::Ptr KdTreePtr
Definition: poisson.h:73
void setMinDepth(int min_depth)
Definition: poisson.h:109
int getIsoDivide()
Get the depth at which a block iso-surface extractor should be used to extract the iso-surface.
Definition: poisson.h:157
shared_ptr< const Poisson< PointNT > > ConstPtr
Definition: poisson.h:65
int getMinDepth()
Definition: poisson.h:112
~Poisson() override
Destructor.
shared_ptr< Poisson< PointNT > > Ptr
Definition: poisson.h:64
void setSolverDivide(int solver_divide)
Set the the depth at which a block Gauss-Seidel solver is used to solve the Laplacian equation.
Definition: poisson.h:140
void setManifold(bool manifold)
Set the manifold flag.
Definition: poisson.h:213
SurfaceReconstruction represents a base surface reconstruction class.
#define PCL_MAKE_ALIGNED_OPERATOR_NEW
Macro to signal a class requires a custom allocator.
Definition: memory.h:63
Defines functions, macros and traits for allocating and using memory.
Defines all the PCL and non-PCL macros used.