pcl::Poisson< PointNT > Class Template Reference
[Module surface]
The Poisson surface reconstruction algorithm.
More...
#include <pcl/surface/poisson.h>
List of all members.
Public Types |
typedef boost::shared_ptr
< Poisson< PointNT > > | Ptr |
typedef boost::shared_ptr
< const Poisson< PointNT > > | ConstPtr |
typedef pcl::PointCloud
< PointNT >::Ptr | PointCloudPtr |
| typedef pcl::KdTree< PointNT > | KdTree |
| typedef pcl::KdTree< PointNT >::Ptr | KdTreePtr |
Public Member Functions |
| | Poisson () |
| | Constructor that sets all the parameters to working default values.
|
| | ~Poisson () |
| | Destructor.
|
| void | performReconstruction (pcl::PolygonMesh &output) |
| | Create the surface.
|
| void | performReconstruction (pcl::PointCloud< PointNT > &points, std::vector< pcl::Vertices > &polygons) |
| | Create the surface.
|
| void | setDepth (int depth) |
| | Set the maximum depth of the tree that will be used for surface reconstruction.
|
| int | getDepth () |
| | Get the depth parameter.
|
| void | setMinDepth (int min_depth) |
| int | getMinDepth () |
| void | setPointWeight (float point_weight) |
| float | getPointWeight () |
| void | setScale (float scale) |
| | Set the ratio between the diameter of the cube used for reconstruction and the diameter of the samples' bounding cube.
|
| float | getScale () |
| | Get the ratio between the diameter of the cube used for reconstruction and the diameter of the samples' bounding cube.
|
| void | setSolverDivide (int solver_divide) |
| | Set the the depth at which a block Gauss-Seidel solver is used to solve the Laplacian equation.
|
| int | getSolverDivide () |
| | Get the the depth at which a block Gauss-Seidel solver is used to solve the Laplacian equation.
|
| void | setIsoDivide (int iso_divide) |
| | Set the depth at which a block iso-surface extractor should be used to extract the iso-surface.
|
| int | getIsoDivide () |
| | Get the depth at which a block iso-surface extractor should be used to extract the iso-surface.
|
| void | setSamplesPerNode (float samples_per_node) |
| | Set the minimum number of sample points that should fall within an octree node as the octree construction is adapted to sampling density.
|
| float | getSamplesPerNode () |
| | Get the minimum number of sample points that should fall within an octree node as the octree construction is adapted to sampling density.
|
| void | setConfidence (bool confidence) |
| | Set the confidence flag.
|
| bool | getConfidence () |
| | Get the confidence flag.
|
| void | setOutputPolygons (bool output_polygons) |
| | Enabling this flag tells the reconstructor to output a polygon mesh (rather than triangulating the results of Marching Cubes).
|
| bool | getOutputPolygons () |
| | Get whether the algorithm outputs a polygon mesh or a triangle mesh.
|
| void | setDegree (int degree) |
| | Set the degree parameter.
|
| int | getDegree () |
| | Get the degree parameter.
|
| void | setManifold (bool manifold) |
| | Set the manifold flag.
|
| bool | getManifold () |
| | Get the manifold flag.
|
Protected Member Functions |
| std::string | getClassName () const |
| | Class get name method.
|
Detailed Description
template<typename PointNT>
class pcl::Poisson< PointNT >
The Poisson surface reconstruction algorithm.
- Note:
- Code adapted from Misha Kazhdan: http://www.cs.jhu.edu/~misha/Code/PoissonRecon/
-
Based on the paper: * Michael Kazhdan, Matthew Bolitho, Hugues Hoppe, "Poisson surface reconstruction", SGP '06 Proceedings of the fourth Eurographics symposium on Geometry processing
- Author:
- Alexandru-Eugen Ichim
Definition at line 60 of file poisson.h.
Member Typedef Documentation
template<typename PointNT >
template<typename PointNT >
template<typename PointNT >
template<typename PointNT >
template<typename PointNT >
Constructor & Destructor Documentation
template<typename PointNT >
Constructor that sets all the parameters to working default values.
Definition at line 64 of file poisson.hpp.
template<typename PointNT >
Member Function Documentation
template<typename PointNT >
| std::string pcl::Poisson< PointNT >::getClassName |
( |
|
) |
const [inline, protected, virtual] |
template<typename PointNT >
Get the confidence flag.
Definition at line 183 of file poisson.h.
template<typename PointNT >
Get the degree parameter.
Definition at line 204 of file poisson.h.
template<typename PointNT >
Get the depth parameter.
Definition at line 105 of file poisson.h.
template<typename PointNT >
Get the depth at which a block iso-surface extractor should be used to extract the iso-surface.
Definition at line 156 of file poisson.h.
template<typename PointNT >
Get the manifold flag.
Definition at line 216 of file poisson.h.
template<typename PointNT >
template<typename PointNT >
| bool pcl::Poisson< PointNT >::getOutputPolygons |
( |
|
) |
[inline] |
Get whether the algorithm outputs a polygon mesh or a triangle mesh.
Definition at line 194 of file poisson.h.
template<typename PointNT >
template<typename PointNT >
| float pcl::Poisson< PointNT >::getSamplesPerNode |
( |
|
) |
[inline] |
Get the minimum number of sample points that should fall within an octree node as the octree construction is adapted to sampling density.
Definition at line 171 of file poisson.h.
template<typename PointNT >
Get the ratio between the diameter of the cube used for reconstruction and the diameter of the samples' bounding cube.
Definition at line 130 of file poisson.h.
template<typename PointNT >
Get the the depth at which a block Gauss-Seidel solver is used to solve the Laplacian equation.
Definition at line 143 of file poisson.h.
template<typename PointNT >
template<typename PointNT >
Create the surface.
- Parameters:
-
| [out] | output | the resultant polygonal mesh |
Implements pcl::SurfaceReconstruction< PointNT >.
Definition at line 152 of file poisson.hpp.
References pcl::PolygonMesh::cloud, pcl::poisson::Point3D< Real >::coords, pcl::poisson::CoredMeshData::inCorePoints, pcl::poisson::CoredVectorMeshData::nextOutOfCorePoint(), pcl::poisson::CoredVectorMeshData::nextPolygon(), pcl::poisson::CoredVectorMeshData::outOfCorePointCount(), pcl::PointCloud< PointT >::points, pcl::poisson::CoredVectorMeshData::polygonCount(), pcl::PolygonMesh::polygons, pcl::toPCLPointCloud2(), and pcl::Vertices::vertices.
template<typename PointNT >
| void pcl::Poisson< PointNT >::setConfidence |
( |
bool |
confidence |
) |
[inline] |
Set the confidence flag.
- Note:
- Enabling this flag tells the reconstructor to use the size of the normals as confidence information. When the flag is not enabled, all normals are normalized to have unit-length prior to reconstruction.
- Parameters:
-
| [in] | confidence | the given flag |
Definition at line 179 of file poisson.h.
template<typename PointNT >
| void pcl::Poisson< PointNT >::setDegree |
( |
int |
degree |
) |
[inline] |
Set the degree parameter.
- Parameters:
-
| [in] | degree | the given degree |
Definition at line 200 of file poisson.h.
template<typename PointNT >
| void pcl::Poisson< PointNT >::setDepth |
( |
int |
depth |
) |
[inline] |
Set the maximum depth of the tree that will be used for surface reconstruction.
- Note:
- Running at depth d corresponds to solving on a voxel grid whose resolution is no larger than 2^d x 2^d x 2^d. Note that since the reconstructor adapts the octree to the sampling density, the specified reconstruction depth is only an upper bound.
- Parameters:
-
| [in] | depth | the depth parameter |
Definition at line 101 of file poisson.h.
template<typename PointNT >
| void pcl::Poisson< PointNT >::setIsoDivide |
( |
int |
iso_divide |
) |
[inline] |
Set the depth at which a block iso-surface extractor should be used to extract the iso-surface.
- Note:
- Using this parameter helps reduce the memory overhead at the cost of a small increase in extraction time. (In practice, we have found that for reconstructions of depth 9 or higher a subdivide depth of 7 or 8 can greatly reduce the memory usage.)
- Parameters:
-
| [in] | iso_divide | the given parameter value |
Definition at line 152 of file poisson.h.
template<typename PointNT >
| void pcl::Poisson< PointNT >::setManifold |
( |
bool |
manifold |
) |
[inline] |
Set the manifold flag.
- Note:
- Enabling this flag tells the reconstructor to add the polygon barycenter when triangulating polygons with more than three vertices.
- Parameters:
-
| [in] | manifold | the given flag |
Definition at line 212 of file poisson.h.
template<typename PointNT >
| void pcl::Poisson< PointNT >::setMinDepth |
( |
int |
min_depth |
) |
[inline] |
template<typename PointNT >
| void pcl::Poisson< PointNT >::setOutputPolygons |
( |
bool |
output_polygons |
) |
[inline] |
Enabling this flag tells the reconstructor to output a polygon mesh (rather than triangulating the results of Marching Cubes).
- Parameters:
-
| [in] | output_polygons | the given flag |
Definition at line 190 of file poisson.h.
template<typename PointNT >
| void pcl::Poisson< PointNT >::setPointWeight |
( |
float |
point_weight |
) |
[inline] |
template<typename PointNT >
| void pcl::Poisson< PointNT >::setSamplesPerNode |
( |
float |
samples_per_node |
) |
[inline] |
Set the minimum number of sample points that should fall within an octree node as the octree construction is adapted to sampling density.
- Note:
- For noise-free samples, small values in the range [1.0 - 5.0] can be used. For more noisy samples, larger values in the range [15.0 - 20.0] may be needed to provide a smoother, noise-reduced, reconstruction.
- Parameters:
-
| [in] | samples_per_node | the given parameter value |
Definition at line 165 of file poisson.h.
template<typename PointNT >
| void pcl::Poisson< PointNT >::setScale |
( |
float |
scale |
) |
[inline] |
Set the ratio between the diameter of the cube used for reconstruction and the diameter of the samples' bounding cube.
- Parameters:
-
| [in] | scale | the given parameter value |
Definition at line 124 of file poisson.h.
template<typename PointNT >
| void pcl::Poisson< PointNT >::setSolverDivide |
( |
int |
solver_divide |
) |
[inline] |
Set the the depth at which a block Gauss-Seidel solver is used to solve the Laplacian equation.
- Note:
- Using this parameter helps reduce the memory overhead at the cost of a small increase in reconstruction time. (In practice, we have found that for reconstructions of depth 9 or higher a subdivide depth of 7 or 8 can greatly reduce the memory usage.)
- Parameters:
-
| [in] | solver_divide | the given parameter value |
Definition at line 139 of file poisson.h.
The documentation for this class was generated from the following files:
- /builddir/build/BUILD/pcl-pcl-1.7.1/surface/include/pcl/surface/poisson.h
- /builddir/build/BUILD/pcl-pcl-1.7.1/surface/include/pcl/surface/impl/poisson.hpp