ATTPCROOT  0.3.0-alpha
A ROOT-based framework for analyzing data from active target detectors
smoothenCloud.cxx
Go to the documentation of this file.
1 #include "smoothenCloud.h"
2 
3 #include <boost/core/checked_delete.hpp> // for checked_delete
4 #include <boost/smart_ptr/shared_ptr.hpp> // for shared_ptr
5 #include <pcl/common/centroid.h>
6 #include <pcl/common/impl/centroid.hpp> // for CentroidPoint::add, Centro...
7 #include <pcl/common/io.h> // for copyPointCloud
8 #include <pcl/kdtree/kdtree_flann.h>
9 
10 #include <stddef.h> // for size_t
11 
12 #include <algorithm> // for nth_element
13 #include <vector> // for vector
14 
15 pcl::PointCloud<pcl::PointXYZI>::Ptr smoothenCloud(pcl::PointCloud<pcl::PointXYZI>::Ptr cloud, double r, bool useMedian)
16 {
17  pcl::PointCloud<pcl::PointXYZI>::Ptr result_cloud(new pcl::PointCloud<pcl::PointXYZI>());
18  if (r == 0) {
19  pcl::copyPointCloud(*cloud, *result_cloud);
20  return result_cloud;
21  }
22  pcl::KdTreeFLANN<pcl::PointXYZI> kdtree;
23 
24  kdtree.setInputCloud(cloud);
25 
26  for (size_t i = 0; i < cloud->size(); ++i) {
27  pcl::PointXYZI newPoint;
28 
29  pcl::CentroidPoint<pcl::PointXYZI> centroid;
30  pcl::PointXYZI centroidPoint;
31  std::vector<int> pointIdxNKNSearch;
32  std::vector<float> pointNKNSquaredDistance;
33  int found = kdtree.radiusSearch(*cloud, (int)i, r, pointIdxNKNSearch, pointNKNSquaredDistance);
34 
35  for (int j = 0; j < found; ++j) {
36  centroid.add((*cloud)[pointIdxNKNSearch[j]]);
37  }
38 
39  centroid.get(centroidPoint);
40 
41  if (useMedian) {
42  std::vector<float> xList;
43  std::vector<float> yList;
44  std::vector<float> zList;
45 
46  for (int j = 0; j < found; ++j) {
47  pcl::PointXYZI const &point = (*cloud)[pointIdxNKNSearch[j]];
48 
49  xList.push_back(point.x);
50  yList.push_back(point.y);
51  zList.push_back(point.z);
52  }
53 
54  auto xMedianIt = xList.begin() + (xList.size() / 2);
55  std::nth_element(xList.begin(), xMedianIt, xList.end());
56  newPoint.x = *xMedianIt;
57 
58  auto yMedianIt = yList.begin() + (yList.size() / 2);
59  std::nth_element(yList.begin(), yMedianIt, yList.end());
60  newPoint.y = *yMedianIt;
61 
62  auto zMedianIt = zList.begin() + (zList.size() / 2);
63  std::nth_element(zList.begin(), zMedianIt, zList.end());
64  newPoint.z = *zMedianIt;
65  } else
66  newPoint = centroidPoint;
67 
68  result_cloud->push_back(newPoint);
69  }
70 
71  return result_cloud;
72 }
smoothenCloud
pcl::PointCloud< pcl::PointXYZI >::Ptr smoothenCloud(pcl::PointCloud< pcl::PointXYZI >::Ptr cloud, double r, bool useMedian)
Definition: smoothenCloud.cxx:15
smoothenCloud.h