Point Cloud Library (PCL)  1.4.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines
concave_hull.hpp
Go to the documentation of this file.
00001 
00040 #include <pcl/pcl_config.h>
00041 #ifdef HAVE_QHULL
00042 
00043 #ifndef PCL_SURFACE_IMPL_CONCAVE_HULL_H_
00044 #define PCL_SURFACE_IMPL_CONCAVE_HULL_H_
00045 
00046 #include <map>
00047 #include <pcl/surface/concave_hull.h>
00048 #include <pcl/common/common.h>
00049 #include <pcl/common/eigen.h>
00050 #include <pcl/common/centroid.h>
00051 #include <pcl/common/transforms.h>
00052 #include <pcl/kdtree/kdtree_flann.h>
00053 #include <pcl/common/io.h>
00054 #include <stdio.h>
00055 #include <stdlib.h>
00056 
00057 extern "C"
00058 {
00059 #ifdef HAVE_QHULL_2011
00060 #  include "libqhull/libqhull.h"
00061 #  include "libqhull/mem.h"
00062 #  include "libqhull/qset.h"
00063 #  include "libqhull/geom.h"
00064 #  include "libqhull/merge.h"
00065 #  include "libqhull/poly.h"
00066 #  include "libqhull/io.h"
00067 #  include "libqhull/stat.h"
00068 #else
00069 #  include "qhull/qhull.h"
00070 #  include "qhull/mem.h"
00071 #  include "qhull/qset.h"
00072 #  include "qhull/geom.h"
00073 #  include "qhull/merge.h"
00074 #  include "qhull/poly.h"
00075 #  include "qhull/io.h"
00076 #  include "qhull/stat.h"
00077 #endif
00078 }
00079 
00081 template <typename PointInT> void
00082 pcl::ConcaveHull<PointInT>::reconstruct (PointCloud &output)
00083 {
00084   output.header = input_->header;
00085   if (alpha_ <= 0)
00086   {
00087     PCL_ERROR ("[pcl::%s::reconstruct] Alpha parameter must be set to a positive number!\n", getClassName ().c_str ());
00088     output.points.clear ();
00089     return;
00090   }
00091 
00092   if (!initCompute ())
00093   {
00094     output.points.clear ();
00095     return;
00096   }
00097 
00098   // Perform the actual surface reconstruction
00099   std::vector<pcl::Vertices> polygons;
00100   performReconstruction (output, polygons);
00101 
00102   output.width = output.points.size ();
00103   output.height = 1;
00104   output.is_dense = true;
00105 
00106   deinitCompute ();
00107 }
00108 
00110 template <typename PointInT> void
00111 pcl::ConcaveHull<PointInT>::reconstruct (PointCloud &output, std::vector<pcl::Vertices> &polygons)
00112 {
00113   output.header = input_->header;
00114   if (alpha_ <= 0)
00115   {
00116     PCL_ERROR ("[pcl::%s::reconstruct] Alpha parameter must be set to a positive number!\n", getClassName ().c_str ());
00117     output.points.clear ();
00118     return;
00119   }
00120 
00121   if (!initCompute ())
00122   {
00123     output.points.clear ();
00124     return;
00125   }
00126 
00127   // Perform the actual surface reconstruction
00128   performReconstruction (output, polygons);
00129 
00130   output.width = output.points.size ();
00131   output.height = 1;
00132   output.is_dense = true;
00133 
00134   deinitCompute ();
00135 }
00136 
00138 template <typename PointInT> void
00139 pcl::ConcaveHull<PointInT>::performReconstruction (PointCloud &alpha_shape, std::vector<pcl::Vertices> &polygons)
00140 {
00141   EIGEN_ALIGN16 Eigen::Matrix3f covariance_matrix;
00142   Eigen::Vector4f xyz_centroid;
00143   compute3DCentroid (*input_, *indices_, xyz_centroid);
00144   computeCovarianceMatrix (*input_, *indices_, xyz_centroid, covariance_matrix);
00145   EIGEN_ALIGN16 Eigen::Vector3f eigen_values;
00146   EIGEN_ALIGN16 Eigen::Matrix3f eigen_vectors;
00147   pcl::eigen33 (covariance_matrix, eigen_vectors, eigen_values);
00148 
00149   Eigen::Affine3f transform1;
00150   transform1.setIdentity ();
00151   int dim = 3;
00152   if (eigen_values[0] / eigen_values[2] < 1.0e-5)
00153   {
00154     // we have points laying on a plane, using 2d convex hull
00155     // compute transformation bring eigen_vectors.col(i) to z-axis
00156 
00157     eigen_vectors.col (2) = eigen_vectors.col (0).cross (eigen_vectors.col (1));
00158     eigen_vectors.col (1) = eigen_vectors.col (2).cross (eigen_vectors.col (0));
00159 
00160     transform1 (0, 2) = eigen_vectors (0, 0);
00161     transform1 (1, 2) = eigen_vectors (1, 0);
00162     transform1 (2, 2) = eigen_vectors (2, 0);
00163 
00164     transform1 (0, 1) = eigen_vectors (0, 1);
00165     transform1 (1, 1) = eigen_vectors (1, 1);
00166     transform1 (2, 1) = eigen_vectors (2, 1);
00167     transform1 (0, 0) = eigen_vectors (0, 2);
00168     transform1 (1, 0) = eigen_vectors (1, 2);
00169     transform1 (2, 0) = eigen_vectors (2, 2);
00170 
00171     transform1 = transform1.inverse ();
00172     dim = 2;
00173   }
00174   else
00175   {
00176     transform1.setIdentity ();
00177   }
00178 
00179   dim_ = dim;
00180 
00181   PointCloud cloud_transformed;
00182   pcl::demeanPointCloud (*input_, *indices_, xyz_centroid, cloud_transformed);
00183   pcl::transformPointCloud (cloud_transformed, cloud_transformed, transform1);
00184 
00185   // True if qhull should free points in qh_freeqhull() or reallocation
00186   boolT ismalloc = True;
00187   // option flags for qhull, see qh_opt.htm
00188   char flags[] = "qhull d QJ";
00189   // output from qh_produce_output(), use NULL to skip qh_produce_output()
00190   FILE *outfile = NULL;
00191   // error messages from qhull code
00192   FILE *errfile = stderr;
00193   // 0 if no error from qhull
00194   int exitcode;
00195 
00196   // Array of coordinates for each point
00197   coordT *points = (coordT*)calloc (cloud_transformed.points.size () * dim, sizeof(coordT));
00198 
00199   for (size_t i = 0; i < cloud_transformed.points.size (); ++i)
00200   {
00201     points[i * dim + 0] = (coordT)cloud_transformed.points[i].x;
00202     points[i * dim + 1] = (coordT)cloud_transformed.points[i].y;
00203 
00204     if (dim > 2)
00205       points[i * dim + 2] = (coordT)cloud_transformed.points[i].z;
00206   }
00207 
00208   // Compute concave hull
00209   exitcode = qh_new_qhull (dim, cloud_transformed.points.size (), points, ismalloc, flags, outfile, errfile);
00210 
00211   if (exitcode != 0)
00212   {
00213     PCL_ERROR ("[pcl::%s::performReconstrution] ERROR: qhull was unable to compute a concave hull for the given point cloud (%lu)!\n", getClassName ().c_str (), (unsigned long) cloud_transformed.points.size ());
00214 
00215     //check if it fails because of NaN values...
00216     if (!cloud_transformed.is_dense)
00217     {
00218       bool NaNvalues = false;
00219       for (size_t i = 0; i < cloud_transformed.size (); ++i)
00220       {
00221         if (!pcl_isfinite (cloud_transformed.points[i].x) || 
00222             !pcl_isfinite (cloud_transformed.points[i].y) ||
00223             !pcl_isfinite (cloud_transformed.points[i].z))
00224         {
00225           NaNvalues = true;
00226           break;
00227         }
00228       }
00229 
00230       if (NaNvalues)
00231         PCL_ERROR ("[pcl::%s::performReconstruction] ERROR: point cloud contains NaN values, consider running pcl::PassThrough filter first to remove NaNs!\n", getClassName ().c_str ());
00232     }
00233 
00234     alpha_shape.points.resize (0);
00235     alpha_shape.width = alpha_shape.height = 0;
00236     polygons.resize (0);
00237 
00238     qh_freeqhull (!qh_ALL);
00239     int curlong, totlong;
00240     qh_memfreeshort (&curlong, &totlong);
00241 
00242     return;
00243   }
00244 
00245   qh_setvoronoi_all ();
00246 
00247   int num_vertices = qh num_vertices;
00248   alpha_shape.points.resize (num_vertices);
00249 
00250   vertexT *vertex;
00251   // Max vertex id
00252   int max_vertex_id = - 1;
00253   FORALLvertices
00254   {
00255     if ((int)vertex->id > max_vertex_id)
00256       max_vertex_id = vertex->id;
00257   }
00258 
00259   facetT *facet;    // set by FORALLfacets
00260 
00261   ++max_vertex_id;
00262   std::vector<int> qhid_to_pcidx (max_vertex_id);
00263 
00264   int num_facets = qh num_facets;
00265   int dd = 0;
00266 
00267   if (dim == 3)
00268   {
00269     setT *triangles_set = qh_settemp (4 * num_facets);
00270     if (voronoi_centers_)
00271       voronoi_centers_->points.resize (num_facets);
00272 
00273     int non_upper = 0;
00274     FORALLfacets
00275     {
00276       // Facets are tetrahedrons (3d)
00277       if (!facet->upperdelaunay)
00278       {
00279         vertexT *anyVertex = (vertexT*)(facet->vertices->e[0].p);
00280         double *center = facet->center;
00281         double r = qh_pointdist (anyVertex->point,center,dim);
00282         facetT *neighb;
00283 
00284         if (voronoi_centers_)
00285         {
00286           voronoi_centers_->points[non_upper].x = facet->center[0];
00287           voronoi_centers_->points[non_upper].y = facet->center[1];
00288           voronoi_centers_->points[non_upper].z = facet->center[2];
00289         }
00290 
00291         non_upper++;
00292 
00293         if (r <= alpha_)
00294         {
00295           // all triangles in tetrahedron are good, add them all to the alpha shape (triangles_set)
00296           qh_makeridges (facet);
00297           facet->good = true;
00298           facet->visitid = qh visit_id;
00299           ridgeT *ridge, **ridgep;
00300           FOREACHridge_ (facet->ridges)
00301           {
00302             neighb = otherfacet_ (ridge, facet);
00303             if ((neighb->visitid != qh visit_id))
00304               qh_setappend (&triangles_set, ridge);
00305           }
00306         }
00307         else
00308         {
00309           // consider individual triangles from the tetrahedron...
00310           facet->good = false;
00311           facet->visitid = qh visit_id;
00312           qh_makeridges (facet);
00313           ridgeT *ridge, **ridgep;
00314           FOREACHridge_ (facet->ridges)
00315           {
00316             facetT *neighb;
00317             neighb = otherfacet_ (ridge, facet);
00318             if ((neighb->visitid != qh visit_id))
00319             {
00320               // check if individual triangle is good and add it to triangles_set
00321 
00322               PointInT a, b, c;
00323               a.x = ((vertexT*)(ridge->vertices->e[0].p))->point[0];
00324               a.y = ((vertexT*)(ridge->vertices->e[0].p))->point[1];
00325               a.z = ((vertexT*)(ridge->vertices->e[0].p))->point[2];
00326               b.x = ((vertexT*)(ridge->vertices->e[1].p))->point[0];
00327               b.y = ((vertexT*)(ridge->vertices->e[1].p))->point[1];
00328               b.z = ((vertexT*)(ridge->vertices->e[1].p))->point[2];
00329               c.x = ((vertexT*)(ridge->vertices->e[2].p))->point[0];
00330               c.y = ((vertexT*)(ridge->vertices->e[2].p))->point[1];
00331               c.z = ((vertexT*)(ridge->vertices->e[2].p))->point[2];
00332 
00333               double r = pcl::getCircumcircleRadius (a, b, c);
00334               if (r <= alpha_)
00335                 qh_setappend (&triangles_set, ridge);
00336             }
00337           }
00338         }
00339       }
00340     }
00341 
00342     if (voronoi_centers_)
00343       voronoi_centers_->points.resize (non_upper);
00344 
00345     // filter, add points to alpha_shape and create polygon structure
00346 
00347     int num_good_triangles = 0;
00348     ridgeT *ridge, **ridgep;
00349     FOREACHridge_ (triangles_set)
00350     {
00351       if (ridge->bottom->upperdelaunay || ridge->top->upperdelaunay || !ridge->top->good || !ridge->bottom->good)
00352         num_good_triangles++;
00353     }
00354 
00355     polygons.resize (num_good_triangles);
00356 
00357     int vertices = 0;
00358     std::vector<bool> added_vertices (max_vertex_id, false);
00359 
00360     int triangles = 0;
00361     FOREACHridge_ (triangles_set)
00362     {
00363       if (ridge->bottom->upperdelaunay || ridge->top->upperdelaunay || !ridge->top->good || !ridge->bottom->good)
00364       {
00365         polygons[triangles].vertices.resize (3);
00366         int vertex_n, vertex_i;
00367         FOREACHvertex_i_ ((*ridge).vertices)  //3 vertices per ridge!
00368         {
00369           if (!added_vertices[vertex->id])
00370           {
00371             alpha_shape.points[vertices].x = vertex->point[0];
00372             alpha_shape.points[vertices].y = vertex->point[1];
00373             alpha_shape.points[vertices].z = vertex->point[2];
00374 
00375             qhid_to_pcidx[vertex->id] = vertices;   //map the vertex id of qhull to the point cloud index
00376             added_vertices[vertex->id] = true;
00377             vertices++;
00378           }
00379 
00380           polygons[triangles].vertices[vertex_i] = qhid_to_pcidx[vertex->id];
00381 
00382         }
00383 
00384         triangles++;
00385       }
00386     }
00387 
00388     alpha_shape.points.resize (vertices);
00389     alpha_shape.width = alpha_shape.points.size ();
00390     alpha_shape.height = 1;
00391   }
00392   else
00393   {
00394     // Compute the alpha complex for the set of points
00395     // Filters the delaunay triangles
00396     setT *edges_set = qh_settemp (3 * num_facets);
00397     if (voronoi_centers_)
00398       voronoi_centers_->points.resize (num_facets);
00399 
00400     FORALLfacets
00401     {
00402       // Facets are the delaunay triangles (2d)
00403       if (!facet->upperdelaunay)
00404       {
00405         // Check if the distance from any vertex to the facet->center
00406         // (center of the voronoi cell) is smaller than alpha
00407         vertexT *anyVertex = (vertexT*)(facet->vertices->e[0].p);
00408         double r = (sqrt ((anyVertex->point[0] - facet->center[0]) *
00409                           (anyVertex->point[0] - facet->center[0]) +
00410                           (anyVertex->point[1] - facet->center[1]) *
00411                           (anyVertex->point[1] - facet->center[1])));
00412         if (r <= alpha_)
00413         {
00414           pcl::Vertices facet_vertices;   //TODO: is not used!!
00415           qh_makeridges (facet);
00416           facet->good = true;
00417 
00418           ridgeT *ridge, **ridgep;
00419           FOREACHridge_ (facet->ridges)
00420           qh_setappend (&edges_set, ridge);
00421 
00422           if (voronoi_centers_)
00423           {
00424             voronoi_centers_->points[dd].x = facet->center[0];
00425             voronoi_centers_->points[dd].y = facet->center[1];
00426             voronoi_centers_->points[dd].z = 0;
00427           }
00428 
00429           ++dd;
00430         }
00431         else
00432           facet->good = false;
00433       }
00434     }
00435 
00436     int vertices = 0;
00437     std::vector<bool> added_vertices (max_vertex_id, false);
00438     std::map<int, std::vector<int> > edges;
00439 
00440     ridgeT *ridge, **ridgep;
00441     FOREACHridge_ (edges_set)
00442     {
00443       if (ridge->bottom->upperdelaunay || ridge->top->upperdelaunay || !ridge->top->good || !ridge->bottom->good)
00444       {
00445         int vertex_n, vertex_i;
00446         int vertices_in_ridge=0;
00447         std::vector<int> pcd_indices;
00448         pcd_indices.resize (2);
00449 
00450         FOREACHvertex_i_ ((*ridge).vertices)  //in 2-dim, 2 vertices per ridge!
00451         {
00452           if (!added_vertices[vertex->id])
00453           {
00454             alpha_shape.points[vertices].x = vertex->point[0];
00455             alpha_shape.points[vertices].y = vertex->point[1];
00456 
00457             if (dim > 2)
00458               alpha_shape.points[vertices].z = vertex->point[2];
00459             else
00460               alpha_shape.points[vertices].z = 0;
00461 
00462             qhid_to_pcidx[vertex->id] = vertices;   //map the vertex id of qhull to the point cloud index
00463             added_vertices[vertex->id] = true;
00464             pcd_indices[vertices_in_ridge] = vertices;
00465             vertices++;
00466           }
00467           else
00468           {
00469             pcd_indices[vertices_in_ridge] = qhid_to_pcidx[vertex->id];
00470           }
00471 
00472           vertices_in_ridge++;
00473         }
00474 
00475         // make edges bidirectional and pointing to alpha_shape pointcloud...
00476         edges[pcd_indices[0]].push_back (pcd_indices[1]);
00477         edges[pcd_indices[1]].push_back (pcd_indices[0]);
00478       }
00479     }
00480 
00481     alpha_shape.points.resize (vertices);
00482 
00483     std::vector<std::vector<int> > connected;
00484     PointCloud alpha_shape_sorted;
00485     alpha_shape_sorted.points.resize (vertices);
00486 
00487     // iterate over edges until they are empty!
00488     std::map<int, std::vector<int> >::iterator curr = edges.begin ();
00489     int next = - 1;
00490     std::vector<bool> used (vertices, false);   // used to decide which direction should we take!
00491     std::vector<int> pcd_idx_start_polygons;
00492     pcd_idx_start_polygons.push_back (0);
00493 
00494     // start following edges and removing elements
00495     int sorted_idx = 0;
00496     while (!edges.empty ())
00497     {
00498       alpha_shape_sorted.points[sorted_idx] = alpha_shape.points[(*curr).first];
00499       // check where we can go from (*curr).first
00500       for (size_t i = 0; i < (*curr).second.size (); i++)
00501       {
00502         if (!used[((*curr).second)[i]])
00503         {
00504           // we can go there
00505           next = ((*curr).second)[i];
00506           break;
00507         }
00508       }
00509 
00510       used[(*curr).first] = true;
00511       edges.erase (curr);   // remove edges starting from curr
00512 
00513       sorted_idx++;
00514 
00515       if (edges.empty ())
00516         break;
00517 
00518       // reassign current
00519       curr = edges.find (next);   // if next is not found, then we have unconnected polygons.
00520       if (curr == edges.end ())
00521       {
00522         // set current to any of the remaining in edge!
00523         curr = edges.begin ();
00524         pcd_idx_start_polygons.push_back (sorted_idx);
00525       }
00526     }
00527 
00528     pcd_idx_start_polygons.push_back (sorted_idx);
00529 
00530     alpha_shape.points = alpha_shape_sorted.points;
00531 
00532     polygons.resize (pcd_idx_start_polygons.size () - 1);
00533 
00534     for (size_t poly_id = 0; poly_id < polygons.size (); poly_id++)
00535     {
00536       polygons[poly_id].vertices.resize (pcd_idx_start_polygons[poly_id + 1] - pcd_idx_start_polygons[poly_id] + 1);
00537       // populate points in the corresponding polygon
00538       for (size_t j = (size_t)pcd_idx_start_polygons[poly_id]; j < (size_t)pcd_idx_start_polygons[poly_id + 1]; ++j)
00539         polygons[poly_id].vertices[j - pcd_idx_start_polygons[poly_id]] = j;
00540 
00541       polygons[poly_id].vertices[polygons[poly_id].vertices.size () - 1] = pcd_idx_start_polygons[poly_id];
00542     }
00543 
00544     if (voronoi_centers_)
00545       voronoi_centers_->points.resize (dd);
00546   }
00547 
00548   qh_freeqhull (!qh_ALL);
00549   int curlong, totlong;
00550   qh_memfreeshort (&curlong, &totlong);
00551 
00552   Eigen::Affine3f transInverse = transform1.inverse ();
00553   pcl::transformPointCloud (alpha_shape, alpha_shape, transInverse);
00554   xyz_centroid[0] = - xyz_centroid[0];
00555   xyz_centroid[1] = - xyz_centroid[1];
00556   xyz_centroid[2] = - xyz_centroid[2];
00557   pcl::demeanPointCloud (alpha_shape, xyz_centroid, alpha_shape);
00558 
00559   // also transform voronoi_centers_...
00560   if (voronoi_centers_)
00561   {
00562     pcl::transformPointCloud (*voronoi_centers_, *voronoi_centers_, transInverse);
00563     pcl::demeanPointCloud (*voronoi_centers_, xyz_centroid, *voronoi_centers_);
00564   }
00565 
00566   if (keep_information_)
00567   {
00568     // build a tree with the original points
00569     pcl::KdTreeFLANN<PointInT> tree (true);
00570     tree.setInputCloud (input_, indices_);
00571 
00572     std::vector<int> neighbor;
00573     std::vector<float> distances;
00574     neighbor.resize (1);
00575     distances.resize (1);
00576 
00577     // for each point in the concave hull, search for the nearest neighbor in the original point cloud
00578 
00579     std::vector<int> indices;
00580     indices.resize (alpha_shape.points.size ());
00581 
00582     for (size_t i = 0; i < alpha_shape.points.size (); i++)
00583     {
00584       tree.nearestKSearch (alpha_shape.points[i], 1, neighbor, distances);
00585       indices[i] = (*indices_)[neighbor[0]];
00586     }
00587 
00588     // replace point with the closest neighbor in the original point cloud
00589     pcl::copyPointCloud (*input_, indices, alpha_shape);
00590   }
00591 }
00592 
00593 template <typename PointInT> void
00594 pcl::ConcaveHull<PointInT>::performReconstruction (PolygonMesh &output)
00595 {  
00596   // Perform reconstruction
00597   pcl::PointCloud<PointInT> hull_points;
00598   performReconstruction (hull_points, output.polygons);
00599 
00600   // Convert the PointCloud into a PointCloud2
00601   pcl::toROSMsg (hull_points, output.cloud);
00602 }
00603 
00604 template <typename PointInT> void
00605 pcl::ConcaveHull<PointInT>::performReconstruction (std::vector<pcl::Vertices> &polygons)
00606 {
00607   pcl::PointCloud<PointInT> hull_points;
00608   performReconstruction (hull_points, polygons);
00609 }
00610 
00611 
00612 #define PCL_INSTANTIATE_ConcaveHull(T) template class PCL_EXPORTS pcl::ConcaveHull<T>;
00613 
00614 #endif    // PCL_SURFACE_IMPL_CONCAVE_HULL_H_
00615 #endif
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines