Point Cloud Library (PCL)  1.13.1-dev
extract_polygonal_prism_data.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 the copyright holder(s) 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 <limits>
41 
42 #include <pcl/pcl_base.h>
43 
44 namespace pcl
45 {
46  /** \brief General purpose method for checking if a 3D point is inside or
47  * outside a given 2D polygon.
48  * \note this method accepts any general 3D point that is projected onto the
49  * 2D polygon, but performs an internal XY projection of both the polygon and the point.
50  * \param point a 3D point projected onto the same plane as the polygon
51  * \param polygon a polygon
52  * \ingroup segmentation
53  */
54  template <typename PointT> bool
55  isPointIn2DPolygon (const PointT &point, const pcl::PointCloud<PointT> &polygon);
56 
57  /** \brief Check if a 2d point (X and Y coordinates considered only!) is
58  * inside or outside a given polygon. This method assumes that both the point
59  * and the polygon are projected onto the XY plane.
60  *
61  * \note (This is highly optimized code taken from http://www.visibone.com/inpoly/)
62  * Copyright (c) 1995-1996 Galacticomm, Inc. Freeware source code.
63  * \param point a 2D point projected onto the same plane as the polygon
64  * \param polygon a polygon
65  * \ingroup segmentation
66  */
67  template <typename PointT> bool
68  isXYPointIn2DXYPolygon (const PointT &point, const pcl::PointCloud<PointT> &polygon);
69 
70  ////////////////////////////////////////////////////////////////////////////////////////////
71  /** \brief @b ExtractPolygonalPrismData uses a set of point indices that
72  * represent a planar model, and together with a given height, generates a 3D
73  * polygonal prism. The polygonal prism is then used to segment all points
74  * lying inside it.
75  *
76  * An example of its usage is to extract the data lying within a set of 3D
77  * boundaries (e.g., objects supported by a plane).
78  *
79  * Example usage:
80  * \code{.cpp}
81  * double z_min = 0., z_max = 0.05; // we want the points above the plane, no farther than 5 cm from the surface
82  * pcl::PointCloud<pcl::PointXYZ>::Ptr hull_points (new pcl::PointCloud<pcl::PointXYZ> ());
83  * pcl::ConvexHull<pcl::PointXYZ> hull;
84  * // hull.setDimension (2); // not necessarily needed, but we need to check the dimensionality of the output
85  * hull.setInputCloud (cloud);
86  * hull.reconstruct (hull_points);
87  * if (hull.getDimension () == 2)
88  * {
89  * pcl::ExtractPolygonalPrismData<pcl::PointXYZ> prism;
90  * prism.setInputCloud (point_cloud);
91  * prism.setInputPlanarHull (hull_points);
92  * prism.setHeightLimits (z_min, z_max);
93  * prism.segment (cloud_indices);
94  * }
95  * else
96  * PCL_ERROR ("The input cloud does not represent a planar surface.\n");
97  * \endcode
98  * \author Radu Bogdan Rusu
99  * \ingroup segmentation
100  */
101  template <typename PointT>
102  class ExtractPolygonalPrismData : public PCLBase<PointT>
103  {
108 
109  public:
111  using PointCloudPtr = typename PointCloud::Ptr;
113 
116 
117  /** \brief Empty constructor. */
119  height_limit_min_ (0),
120  height_limit_max_(std::numeric_limits<float>::max()),
121  vpx_ (0), vpy_ (0), vpz_ (0)
122  {};
123 
124  /** \brief Provide a pointer to the input planar hull dataset.
125  * \note Please see the example in the class description for how to obtain this.
126  * \param[in] hull the input planar hull dataset
127  */
128  inline void
130 
131  /** \brief Get a pointer the input planar hull dataset. */
132  inline PointCloudConstPtr
133  getInputPlanarHull () const { return (planar_hull_); }
134 
135  /** \brief Set the height limits. All points having distances to the
136  * model outside this interval will be discarded.
137  *
138  * \param[in] height_min the minimum allowed distance to the plane model value
139  * \param[in] height_max the maximum allowed distance to the plane model value
140  */
141  inline void
142  setHeightLimits (double height_min, double height_max)
143  {
144  height_limit_min_ = height_min;
145  height_limit_max_ = height_max;
146  }
147 
148  /** \brief Get the height limits (min/max) as set by the user. The
149  * default values are 0, std::numeric_limits<float>::max().
150  * \param[out] height_min the resultant min height limit
151  * \param[out] height_max the resultant max height limit
152  */
153  inline void
154  getHeightLimits (double &height_min, double &height_max) const
155  {
156  height_min = height_limit_min_;
157  height_max = height_limit_max_;
158  }
159 
160  /** \brief Set the viewpoint.
161  * \param[in] vpx the X coordinate of the viewpoint
162  * \param[in] vpy the Y coordinate of the viewpoint
163  * \param[in] vpz the Z coordinate of the viewpoint
164  */
165  inline void
166  setViewPoint (float vpx, float vpy, float vpz)
167  {
168  vpx_ = vpx;
169  vpy_ = vpy;
170  vpz_ = vpz;
171  }
172 
173  /** \brief Get the viewpoint. */
174  inline void
175  getViewPoint (float &vpx, float &vpy, float &vpz) const
176  {
177  vpx = vpx_;
178  vpy = vpy_;
179  vpz = vpz_;
180  }
181 
182  /** \brief Cluster extraction in a PointCloud given by <setInputCloud (), setIndices ()>
183  * \param[out] output the resultant point indices that support the model found (inliers)
184  */
185  void
186  segment (PointIndices &output);
187 
188  protected:
189  /** \brief A pointer to the input planar hull dataset. */
191 
192  /** \brief The minimum number of points needed on the convex hull. */
194 
195  /** \brief The minimum allowed height (distance to the model) a point
196  * will be considered from.
197  */
199 
200  /** \brief The maximum allowed height (distance to the model) a point
201  * will be considered from.
202  */
204 
205  /** \brief Values describing the data acquisition viewpoint. Default: 0,0,0. */
206  float vpx_, vpy_, vpz_;
207 
208  /** \brief Class getName method. */
209  virtual std::string
210  getClassName () const { return ("ExtractPolygonalPrismData"); }
211  };
212 }
213 
214 #ifdef PCL_NO_PRECOMPILE
215 #include <pcl/segmentation/impl/extract_polygonal_prism_data.hpp>
216 #endif
ExtractPolygonalPrismData uses a set of point indices that represent a planar model,...
PointCloudConstPtr getInputPlanarHull() const
Get a pointer the input planar hull dataset.
void setHeightLimits(double height_min, double height_max)
Set the height limits.
void setInputPlanarHull(const PointCloudConstPtr &hull)
Provide a pointer to the input planar hull dataset.
void getHeightLimits(double &height_min, double &height_max) const
Get the height limits (min/max) as set by the user.
void getViewPoint(float &vpx, float &vpy, float &vpz) const
Get the viewpoint.
int min_pts_hull_
The minimum number of points needed on the convex hull.
double height_limit_min_
The minimum allowed height (distance to the model) a point will be considered from.
virtual std::string getClassName() const
Class getName method.
float vpx_
Values describing the data acquisition viewpoint.
void setViewPoint(float vpx, float vpy, float vpz)
Set the viewpoint.
double height_limit_max_
The maximum allowed height (distance to the model) a point will be considered from.
PointCloudConstPtr planar_hull_
A pointer to the input planar hull dataset.
void segment(PointIndices &output)
Cluster extraction in a PointCloud given by <setInputCloud (), setIndices ()>
PCL base class.
Definition: pcl_base.h:70
typename PointCloud::Ptr PointCloudPtr
Definition: pcl_base.h:73
typename PointCloud::ConstPtr PointCloudConstPtr
Definition: pcl_base.h:74
PointIndices::ConstPtr PointIndicesConstPtr
Definition: pcl_base.h:77
PointIndices::Ptr PointIndicesPtr
Definition: pcl_base.h:76
PointCloud represents the base class in PCL for storing collections of 3D points.
Definition: point_cloud.h:173
shared_ptr< PointCloud< PointT > > Ptr
Definition: point_cloud.h:413
shared_ptr< const PointCloud< PointT > > ConstPtr
Definition: point_cloud.h:414
bool isPointIn2DPolygon(const PointT &point, const pcl::PointCloud< PointT > &polygon)
General purpose method for checking if a 3D point is inside or outside a given 2D polygon.
bool isXYPointIn2DXYPolygon(const PointT &point, const pcl::PointCloud< PointT > &polygon)
Check if a 2d point (X and Y coordinates considered only!) is inside or outside a given polygon.
shared_ptr< ::pcl::PointIndices > Ptr
Definition: PointIndices.h:13
shared_ptr< const ::pcl::PointIndices > ConstPtr
Definition: PointIndices.h:14
A point structure representing Euclidean xyz coordinates, and the RGB color.