Point Cloud Library (PCL)  1.4.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines
convex_hull.hpp
Go to the documentation of this file.
00001 /*
00002  * Software License Agreement (BSD License)
00003  *
00004  *  Copyright (c) 2010, Willow Garage, Inc.
00005  *  All rights reserved.
00006  *
00007  *  Redistribution and use in source and binary forms, with or without
00008  *  modification, are permitted provided that the following conditions
00009  *  are met:
00010  *
00011  *   * Redistributions of source code must retain the above copyright
00012  *     notice, this list of conditions and the following disclaimer.
00013  *   * Redistributions in binary form must reproduce the above
00014  *     copyright notice, this list of conditions and the following
00015  *     disclaimer in the documentation and/or other materials provided
00016  *     with the distribution.
00017  *   * Neither the name of Willow Garage, Inc. nor the names of its
00018  *     contributors may be used to endorse or promote products derived
00019  *     from this software without specific prior written permission.
00020  *
00021  *  THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
00022  *  "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
00023  *  LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
00024  *  FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
00025  *  COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
00026  *  INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
00027  *  BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
00028  *  LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
00029  *  CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
00030  *  LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
00031  *  ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
00032  *  POSSIBILITY OF SUCH DAMAGE.
00033  *
00034  * $Id: convex_hull.hpp 3753 2011-12-31 23:30:57Z rusu $
00035  *
00036  */
00037 
00038 #include <pcl/pcl_config.h>
00039 #ifdef HAVE_QHULL
00040 
00041 #ifndef PCL_SURFACE_IMPL_CONVEX_HULL_H_
00042 #define PCL_SURFACE_IMPL_CONVEX_HULL_H_
00043 
00044 #include "pcl/surface/convex_hull.h"
00045 #include <pcl/common/common.h>
00046 #include <pcl/common/eigen.h>
00047 #include <pcl/common/transforms.h>
00048 #include <pcl/common/io.h>
00049 #include <pcl/kdtree/kdtree.h>
00050 #include <pcl/kdtree/kdtree_flann.h>
00051 #include <stdio.h>
00052 #include <stdlib.h>
00053 
00054 extern "C"
00055 {
00056 #ifdef HAVE_QHULL_2011
00057 #  include "libqhull/libqhull.h"
00058 #  include "libqhull/mem.h"
00059 #  include "libqhull/qset.h"
00060 #  include "libqhull/geom.h"
00061 #  include "libqhull/merge.h"
00062 #  include "libqhull/poly.h"
00063 #  include "libqhull/io.h"
00064 #  include "libqhull/stat.h"
00065 #else
00066 #  include "qhull/qhull.h"
00067 #  include "qhull/mem.h"
00068 #  include "qhull/qset.h"
00069 #  include "qhull/geom.h"
00070 #  include "qhull/merge.h"
00071 #  include "qhull/poly.h"
00072 #  include "qhull/io.h"
00073 #  include "qhull/stat.h"
00074 #endif
00075 }
00076 
00078 template <typename PointInT> void
00079 pcl::ConvexHull<PointInT>::performReconstruction (PointCloud &hull, std::vector<pcl::Vertices> &polygons,
00080                                                   bool fill_polygon_data)
00081 {
00082   // FInd the principal directions
00083   EIGEN_ALIGN16 Eigen::Matrix3f covariance_matrix;
00084   Eigen::Vector4f xyz_centroid;
00085   compute3DCentroid (*input_, *indices_, xyz_centroid);
00086   computeCovarianceMatrix (*input_, *indices_, xyz_centroid, covariance_matrix);
00087   EIGEN_ALIGN16 Eigen::Vector3f eigen_values;
00088   EIGEN_ALIGN16 Eigen::Matrix3f eigen_vectors;
00089   pcl::eigen33 (covariance_matrix, eigen_vectors, eigen_values);
00090 
00091   Eigen::Affine3f transform1;
00092   transform1.setIdentity ();
00093   int dim = 3;
00094 
00095   if (eigen_values[0] / eigen_values[2] < 1.0e-5)
00096   {
00097     // We have points laying on a plane, using 2d convex hull
00098     // Compute transformation bring eigen_vectors.col(i) to z-axis
00099 
00100     eigen_vectors.col (2) = eigen_vectors.col (0).cross (eigen_vectors.col (1));
00101     eigen_vectors.col (1) = eigen_vectors.col (2).cross (eigen_vectors.col (0));
00102 
00103     transform1 (0, 2) = eigen_vectors (0, 0);
00104     transform1 (1, 2) = eigen_vectors (1, 0);
00105     transform1 (2, 2) = eigen_vectors (2, 0);
00106 
00107     transform1 (0, 1) = eigen_vectors (0, 1);
00108     transform1 (1, 1) = eigen_vectors (1, 1);
00109     transform1 (2, 1) = eigen_vectors (2, 1);
00110     transform1 (0, 0) = eigen_vectors (0, 2);
00111     transform1 (1, 0) = eigen_vectors (1, 2);
00112     transform1 (2, 0) = eigen_vectors (2, 2);
00113 
00114     transform1 = transform1.inverse ();
00115     dim = 2;
00116   }
00117   else
00118     transform1.setIdentity ();
00119 
00120   dim_ = dim;
00121 
00122   PointCloud cloud_transformed;
00123   pcl::demeanPointCloud (*input_, *indices_, xyz_centroid, cloud_transformed);
00124   pcl::transformPointCloud (cloud_transformed, cloud_transformed, transform1);
00125 
00126   // True if qhull should free points in qh_freeqhull() or reallocation
00127   boolT ismalloc = True;
00128   // output from qh_produce_output(), use NULL to skip qh_produce_output()
00129   FILE *outfile = NULL;
00130 
00131   std::string flags_str;
00132   flags_str = "qhull Tc";
00133 
00134   if (compute_area_)
00135   {
00136     flags_str.append (" FA");
00137     outfile = stderr;
00138   }
00139 
00140   // option flags for qhull, see qh_opt.htm
00141   //char flags[] = "qhull Tc FA";
00142   char * flags = (char *)flags_str.c_str();
00143   // error messages from qhull code
00144   FILE *errfile = stderr;
00145   // 0 if no error from qhull
00146   int exitcode;
00147 
00148   // Array of coordinates for each point
00149   coordT *points = (coordT *)calloc (cloud_transformed.points.size () * dim, sizeof(coordT));
00150 
00151   for (size_t i = 0; i < cloud_transformed.points.size (); ++i)
00152   {
00153     points[i * dim + 0] = (coordT)cloud_transformed.points[i].x;
00154     points[i * dim + 1] = (coordT)cloud_transformed.points[i].y;
00155 
00156     if (dim > 2)
00157       points[i * dim + 2] = (coordT)cloud_transformed.points[i].z;
00158   }
00159 
00160   // Compute convex hull
00161   exitcode = qh_new_qhull (dim, cloud_transformed.points.size (), points, ismalloc, flags, outfile, errfile);
00162 
00163   if (exitcode != 0)
00164   {
00165     PCL_ERROR ("[pcl::%s::performReconstrution] ERROR: qhull was unable to compute a convex hull for the given point cloud (%lu)!\n", getClassName ().c_str (), (unsigned long) input_->points.size ());
00166 
00167     //check if it fails because of NaN values...
00168     if (!cloud_transformed.is_dense)
00169     {
00170       bool NaNvalues = false;
00171       for (size_t i = 0; i < cloud_transformed.size (); ++i)
00172       {
00173         if (!pcl_isfinite (cloud_transformed.points[i].x) || 
00174             !pcl_isfinite (cloud_transformed.points[i].y) ||
00175             !pcl_isfinite (cloud_transformed.points[i].z))
00176         {
00177           NaNvalues = true;
00178           break;
00179         }
00180       }
00181 
00182       if (NaNvalues)
00183         PCL_ERROR ("[pcl::%s::performReconstruction] ERROR: point cloud contains NaN values, consider running pcl::PassThrough filter first to remove NaNs!\n", getClassName ().c_str ());
00184     }
00185 
00186     hull.points.resize (0);
00187     hull.width = hull.height = 0;
00188     polygons.resize (0);
00189 
00190     qh_freeqhull (!qh_ALL);
00191     int curlong, totlong;
00192     qh_memfreeshort (&curlong, &totlong);
00193 
00194     return;
00195   }
00196 
00197   qh_triangulate ();
00198 
00199   int num_facets = qh num_facets;
00200 
00201   int num_vertices = qh num_vertices;
00202   hull.points.resize (num_vertices);
00203 
00204   vertexT * vertex;
00205   int i = 0;
00206   // Max vertex id
00207   int max_vertex_id = -1;
00208   FORALLvertices
00209   {
00210     if ((int)vertex->id > max_vertex_id)
00211     max_vertex_id = vertex->id;
00212   }
00213 
00214   ++max_vertex_id;
00215   std::vector<int> qhid_to_pcidx (max_vertex_id);
00216 
00217   FORALLvertices
00218   {
00219     // Add vertices to hull point_cloud
00220     hull.points[i].x = vertex->point[0];
00221     hull.points[i].y = vertex->point[1];
00222 
00223     if (dim>2)
00224     hull.points[i].z = vertex->point[2];
00225     else
00226     hull.points[i].z = 0;
00227 
00228     qhid_to_pcidx[vertex->id] = i; // map the vertex id of qhull to the point cloud index
00229     ++i;
00230   }
00231 
00232   if (compute_area_)
00233   {
00234     total_area_  = qh totarea;
00235     total_volume_ = qh totvol;
00236   }
00237 
00238   if (fill_polygon_data)
00239   {
00240     if (dim == 3)
00241     {
00242       polygons.resize (num_facets);
00243       int dd = 0;
00244 
00245       facetT * facet;
00246       FORALLfacets
00247       {
00248         polygons[dd].vertices.resize (3);
00249 
00250         // Needed by FOREACHvertex_i_
00251         int vertex_n, vertex_i;
00252         FOREACHvertex_i_((*facet).vertices)
00253         //facet_vertices.vertices.push_back (qhid_to_pcidx[vertex->id]);
00254         polygons[dd].vertices[vertex_i] = qhid_to_pcidx[vertex->id];
00255         ++dd;
00256       }
00257 
00258 
00259     }
00260     else
00261     {
00262       // dim=2, we want to return just a polygon with all vertices sorted
00263       // so that they form a non-intersecting polygon...
00264       // this can be moved to the upper loop probably and here just
00265       // the sorting + populate
00266 
00267       Eigen::Vector4f centroid;
00268       pcl::compute3DCentroid (hull, centroid);
00269       centroid[3] = 0;
00270       polygons.resize (1);
00271 
00272       int num_vertices = qh num_vertices, dd = 0;
00273 
00274       // Copy all vertices
00275       std::vector<std::pair<int, Eigen::Vector4f>, Eigen::aligned_allocator<std::pair<int, Eigen::Vector4f> > > idx_points (num_vertices);
00276 
00277       FORALLvertices
00278       {
00279         idx_points[dd].first = qhid_to_pcidx[vertex->id];
00280         idx_points[dd].second = hull.points[idx_points[dd].first].getVector4fMap () - centroid;
00281         ++dd;
00282       }
00283 
00284       // Sort idx_points
00285       std::sort (idx_points.begin (), idx_points.end (), comparePoints2D);
00286       polygons[0].vertices.resize (idx_points.size () + 1);
00287 
00288       //Sort also points...
00289       PointCloud hull_sorted;
00290       hull_sorted.points.resize (hull.points.size ());
00291 
00292       for (size_t j = 0; j < idx_points.size (); ++j)
00293       hull_sorted.points[j] = hull.points[idx_points[j].first];
00294 
00295       hull.points = hull_sorted.points;
00296 
00297       // Populate points
00298       for (size_t j = 0; j < idx_points.size (); ++j)
00299       polygons[0].vertices[j] = j;
00300 
00301       polygons[0].vertices[idx_points.size ()] = 0;
00302     }
00303   }
00304   else
00305   {
00306     if (dim == 2)
00307     {
00308       // We want to sort the points
00309       Eigen::Vector4f centroid;
00310       pcl::compute3DCentroid (hull, centroid);
00311       centroid[3] = 0;
00312       polygons.resize (1);
00313 
00314       int num_vertices = qh num_vertices, dd = 0;
00315 
00316       // Copy all vertices
00317       std::vector<std::pair<int, Eigen::Vector4f>, Eigen::aligned_allocator<std::pair<int, Eigen::Vector4f> > > idx_points (num_vertices);
00318 
00319       FORALLvertices
00320       {
00321         idx_points[dd].first = qhid_to_pcidx[vertex->id];
00322         idx_points[dd].second = hull.points[idx_points[dd].first].getVector4fMap () - centroid;
00323         ++dd;
00324       }
00325 
00326       // Sort idx_points
00327       std::sort (idx_points.begin (), idx_points.end (), comparePoints2D);
00328 
00329       //Sort also points...
00330       PointCloud hull_sorted;
00331       hull_sorted.points.resize (hull.points.size ());
00332 
00333       for (size_t j = 0; j < idx_points.size (); ++j)
00334       hull_sorted.points[j] = hull.points[idx_points[j].first];
00335 
00336       hull.points = hull_sorted.points;
00337     }
00338   }
00339 
00340   // Deallocates memory (also the points)
00341   qh_freeqhull(!qh_ALL);
00342   int curlong, totlong;
00343   qh_memfreeshort (&curlong, &totlong);
00344 
00345   // Rotate the hull point cloud by transform's inverse
00346   // If the input point cloud has been rotated
00347   if (dim == 2)
00348   {
00349     Eigen::Affine3f transInverse = transform1.inverse ();
00350     pcl::transformPointCloud (hull, hull, transInverse);
00351 
00352     //for 2D sets, the qhull library delivers the actual area of the 2d hull in the volume
00353     if(compute_area_)
00354     {
00355       total_area_ = total_volume_;
00356       total_volume_ = 0.0;
00357     }
00358 
00359   }
00360 
00361   xyz_centroid[0] = -xyz_centroid[0];
00362   xyz_centroid[1] = -xyz_centroid[1];
00363   xyz_centroid[2] = -xyz_centroid[2];
00364   pcl::demeanPointCloud (hull, xyz_centroid, hull);
00365 
00366   //if keep_information_
00367   if (keep_information_)
00368   {
00369     //build a tree with the original points
00370     pcl::KdTreeFLANN<PointInT> tree (true);
00371     tree.setInputCloud (input_, indices_);
00372 
00373     std::vector<int> neighbor;
00374     std::vector<float> distances;
00375     neighbor.resize (1);
00376     distances.resize (1);
00377 
00378     //for each point in the convex hull, search for the nearest neighbor
00379 
00380     std::vector<int> indices;
00381     indices.resize (hull.points.size ());
00382 
00383     for (size_t i = 0; i < hull.points.size (); i++)
00384     {
00385       tree.nearestKSearch (hull.points[i], 1, neighbor, distances);
00386       indices[i] = (*indices_)[neighbor[0]];
00387     }
00388 
00389     //replace point with the closest neighbor in the original point cloud
00390     pcl::copyPointCloud (*input_, indices, hull);
00391   }
00392 
00393   hull.width = hull.points.size ();
00394   hull.height = 1;
00395   hull.is_dense = false;
00396 }
00397 
00399 template <typename PointInT> void
00400 pcl::ConvexHull<PointInT>::reconstruct (PointCloud &output)
00401 {
00402   output.header = input_->header;
00403   if (!initCompute () || input_->points.empty ())
00404   {
00405     output.points.clear ();
00406     return;
00407   }
00408 
00409   // Perform the actual surface reconstruction
00410   std::vector<pcl::Vertices> polygons;
00411   performReconstruction (output, polygons, false);
00412 
00413   output.width = output.points.size ();
00414   output.height = 1;
00415   output.is_dense = true;
00416 
00417   deinitCompute ();
00418 }
00419 
00420 
00421 template <typename PointInT> void
00422 pcl::ConvexHull<PointInT>::performReconstruction (PolygonMesh &output)
00423 {  
00424   // Perform reconstruction
00425   pcl::PointCloud<PointInT> hull_points;
00426   performReconstruction (hull_points, output.polygons, true);
00427 
00428   // Convert the PointCloud into a PointCloud2
00429   pcl::toROSMsg (hull_points, output.cloud);
00430 }
00431 
00432 template <typename PointInT> void
00433 pcl::ConvexHull<PointInT>::performReconstruction (std::vector<pcl::Vertices> &polygons)
00434 {
00435   pcl::PointCloud<PointInT> hull_points;
00436   performReconstruction (hull_points, polygons, true);
00437 }
00438 
00440 template <typename PointInT> void
00441 pcl::ConvexHull<PointInT>::reconstruct (PointCloud &points, std::vector<pcl::Vertices> &polygons)
00442 {
00443   points.header = input_->header;
00444   if (!initCompute () || input_->points.empty ())
00445   {
00446     points.points.clear ();
00447     return;
00448   }
00449 
00450   // Perform the actual surface reconstruction
00451   performReconstruction (points, polygons, true);
00452 
00453   points.width = points.points.size ();
00454   points.height = 1;
00455   points.is_dense = true;
00456 
00457   deinitCompute ();
00458 }
00459 
00460 #define PCL_INSTANTIATE_ConvexHull(T) template class PCL_EXPORTS pcl::ConvexHull<T>;
00461 
00462 #endif    // PCL_SURFACE_IMPL_CONVEX_HULL_H_
00463 #endif
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines