Point Cloud Library (PCL) 1.15.1-dev
Loading...
Searching...
No Matches
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
44namespace 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
67 using SurfaceReconstruction<PointNT>::input_;
68 using SurfaceReconstruction<PointNT>::tree_;
69
71
72 /** \brief Constructor that sets all the parameters to working default values. */
74
75 /** \brief Destructor. */
76 ~Poisson () override;
77
78 /** \brief Create the surface.
79 * \param[out] output the resultant polygonal mesh
80 */
81 void
82 performReconstruction (pcl::PolygonMesh &output) override;
83
84 /** \brief Create the surface.
85 * \param[out] points the vertex positions of the resulting mesh
86 * \param[out] polygons the connectivity of the resulting mesh
87 */
88 void
90 std::vector<pcl::Vertices> &polygons) override;
91
92 /** \brief Set the maximum depth of the tree that will be used for surface reconstruction.
93 * \note Running at depth d corresponds to solving on a voxel grid whose resolution is no larger than
94 * 2^d x 2^d x 2^d. Note that since the reconstructor adapts the octree to the sampling density, the specified
95 * reconstruction depth is only an upper bound.
96 * \param[in] depth the depth parameter
97 */
98 inline void
99 setDepth (int depth) { depth_ = depth; }
100
101 /** \brief Get the depth parameter */
102 inline int
103 getDepth () { return depth_; }
104
105 inline void
106 setMinDepth (int min_depth) { min_depth_ = min_depth; }
107
108 inline int
109 getMinDepth () { return min_depth_; }
110
111 inline void
112 setPointWeight (float point_weight) { point_weight_ = point_weight; }
113
114 inline float
115 getPointWeight () { return point_weight_; }
116
117 /** \brief Set the ratio between the diameter of the cube used for reconstruction and the diameter of the
118 * samples' bounding cube.
119 * \param[in] scale the given parameter value
120 */
121 inline void
122 setScale (float scale) { scale_ = scale; }
123
124 /** Get the ratio between the diameter of the cube used for reconstruction and the diameter of the
125 * samples' bounding cube.
126 */
127 inline float
128 getScale () { return scale_; }
129
130 /** \brief Set the the depth at which a block Gauss-Seidel solver is used to solve the Laplacian equation
131 * \note Using this parameter helps reduce the memory overhead at the cost of a small increase in
132 * reconstruction time. (In practice, we have found that for reconstructions of depth 9 or higher a subdivide
133 * depth of 7 or 8 can greatly reduce the memory usage.)
134 * \param[in] solver_divide the given parameter value
135 */
136 inline void
137 setSolverDivide (int solver_divide) { solver_divide_ = solver_divide; }
138
139 /** \brief Get the depth at which a block Gauss-Seidel solver is used to solve the Laplacian equation */
140 inline int
141 getSolverDivide () { return solver_divide_; }
142
143 /** \brief Set the depth at which a block iso-surface extractor should be used to extract the iso-surface
144 * \note Using this parameter helps reduce the memory overhead at the cost of a small increase in extraction
145 * time. (In practice, we have found that for reconstructions of depth 9 or higher a subdivide depth of 7 or 8
146 * can greatly reduce the memory usage.)
147 * \param[in] iso_divide the given parameter value
148 */
149 inline void
150 setIsoDivide (int iso_divide) { iso_divide_ = iso_divide; }
151
152 /** \brief Get the depth at which a block iso-surface extractor should be used to extract the iso-surface */
153 inline int
154 getIsoDivide () { return iso_divide_; }
155
156 /** \brief Set the minimum number of sample points that should fall within an octree node as the octree
157 * construction is adapted to sampling density
158 * \note For noise-free samples, small values in the range [1.0 - 5.0] can be used. For more noisy samples,
159 * larger values in the range [15.0 - 20.0] may be needed to provide a smoother, noise-reduced, reconstruction.
160 * \param[in] samples_per_node the given parameter value
161 */
162 inline void
163 setSamplesPerNode (float samples_per_node) { samples_per_node_ = samples_per_node; }
164
165 /** \brief Get the minimum number of sample points that should fall within an octree node as the octree
166 * construction is adapted to sampling density
167 */
168 inline float
169 getSamplesPerNode () { return samples_per_node_; }
170
171 /** \brief Set the confidence flag
172 * \note Enabling this flag tells the reconstructor to use the size of the normals as confidence information.
173 * When the flag is not enabled, all normals are normalized to have unit-length prior to reconstruction.
174 * \param[in] confidence the given flag
175 */
176 inline void
177 setConfidence (bool confidence) { confidence_ = confidence; }
178
179 /** \brief Get the confidence flag */
180 inline bool
181 getConfidence () { return confidence_; }
182
183 /** \brief Enabling this flag tells the reconstructor to output a polygon mesh (rather than triangulating the
184 * results of Marching Cubes).
185 * \param[in] output_polygons the given flag
186 */
187 inline void
188 setOutputPolygons (bool output_polygons) { output_polygons_ = output_polygons; }
189
190 /** \brief Get whether the algorithm outputs a polygon mesh or a triangle mesh */
191 inline bool
192 getOutputPolygons () { return output_polygons_; }
193
194 /** \brief Set the degree parameter
195 * \param[in] degree the given degree
196 */
197 inline void
198 setDegree (int degree) { degree_ = degree; }
199
200 /** \brief Get the degree parameter */
201 inline int
202 getDegree () { return degree_; }
203
204 /** \brief Set the manifold flag.
205 * \note Enabling this flag tells the reconstructor to add the polygon barycenter when triangulating polygons
206 * with more than three vertices.
207 * \param[in] manifold the given flag
208 */
209 inline void
210 setManifold (bool manifold) { manifold_ = manifold; }
211
212 /** \brief Get the manifold flag */
213 inline bool
214 getManifold () { return manifold_; }
215
216 /** \brief Set the number of threads to use.
217 * \param[in] threads the number of threads
218 */
219 void
220 setThreads(int threads);
221
222
223 /** \brief Get the number of threads*/
224 inline int
226 {
227 return threads_;
228 }
229
230 protected:
231 /** \brief Class get name method. */
232 std::string
233 getClassName () const override { return ("Poisson"); }
234
235 private:
236 int depth_{8};
237 int min_depth_{5};
238 float point_weight_{4};
239 float scale_{1.1f};
240 int solver_divide_{8};
241 int iso_divide_{8};
242 float samples_per_node_{1.0};
243 bool confidence_{false};
244 bool output_polygons_{false};
245
246 bool no_reset_samples_{false};
247 bool no_clip_tree_{false};
248 bool manifold_{true};
249
250 int refine_{3};
251 int kernel_depth_{8};
252 int degree_{2};
253 bool non_adaptive_weights_{false};
254 bool show_residual_{false};
255 int min_iterations_{8};
256 float solver_accuracy_{1e-3f};
257 int threads_{1};
258
259 template<int Degree> void
260 execute (poisson::CoredVectorMeshData &mesh,
261 poisson::Point3D<float> &translate,
262 float &scale);
263
264 public:
266 };
267}
268
269#ifdef PCL_NO_PRECOMPILE
270#include <pcl/surface/impl/poisson.hpp>
271#endif
PointCloudConstPtr input_
The input point cloud dataset.
Definition pcl_base.h:147
KdTreePtr tree_
A pointer to the spatial search object.
PointCloud represents the base class in PCL for storing collections of 3D points.
shared_ptr< PointCloud< PointT > > Ptr
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:233
int getThreads()
Get the number of threads.
Definition poisson.h:225
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:163
float getPointWeight()
Definition poisson.h:115
void setDegree(int degree)
Set the degree parameter.
Definition poisson.h:198
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:103
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:150
float getSamplesPerNode()
Get the minimum number of sample points that should fall within an octree node as the octree construc...
Definition poisson.h:169
void setConfidence(bool confidence)
Set the confidence flag.
Definition poisson.h:177
bool getConfidence()
Get the confidence flag.
Definition poisson.h:181
void setPointWeight(float point_weight)
Definition poisson.h:112
int getSolverDivide()
Get the depth at which a block Gauss-Seidel solver is used to solve the Laplacian equation.
Definition poisson.h:141
int getDegree()
Get the degree parameter.
Definition poisson.h:202
void setOutputPolygons(bool output_polygons)
Enabling this flag tells the reconstructor to output a polygon mesh (rather than triangulating the re...
Definition poisson.h:188
void setDepth(int depth)
Set the maximum depth of the tree that will be used for surface reconstruction.
Definition poisson.h:99
bool getOutputPolygons()
Get whether the algorithm outputs a polygon mesh or a triangle mesh.
Definition poisson.h:192
float getScale()
Get the ratio between the diameter of the cube used for reconstruction and the diameter of the sample...
Definition poisson.h:128
bool getManifold()
Get the manifold flag.
Definition poisson.h:214
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:122
void setMinDepth(int min_depth)
Definition poisson.h:106
int getIsoDivide()
Get the depth at which a block iso-surface extractor should be used to extract the iso-surface.
Definition poisson.h:154
shared_ptr< const Poisson< PointNT > > ConstPtr
Definition poisson.h:65
int getMinDepth()
Definition poisson.h:109
~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:137
void setManifold(bool manifold)
Set the manifold flag.
Definition poisson.h:210
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:86
Defines functions, macros and traits for allocating and using memory.
Defines all the PCL and non-PCL macros used.