ATTPCROOT  0.3.0-alpha
A ROOT-based framework for analyzing data from active target detectors
hc.h
Go to the documentation of this file.
1 #ifndef HC_H
2 #define HC_H
3 
4 #include <Eigen/Core> // for MatrixXf, Vector3f, DenseCoeffsBase
5 #include <pcl/point_cloud.h> // for PointCloud, PointCloud<>::ConstPtr
6 #include <pcl/point_types.h> // for PointXYZI
7 
8 #include "cluster.h"
9 
10 #include <cstddef> // for size_t
11 #include <limits>
12 #include <vector>
13 
14 namespace hc {
15 struct triplet {
16  size_t pointIndexA;
17  size_t pointIndexB;
18  size_t pointIndexC;
19  Eigen::Vector3f center;
20  Eigen::Vector3f direction;
21  float error;
22  friend bool operator<(const triplet &t1, const triplet &t2) { return (t1.error < t2.error); };
23 };
24 
25 typedef std::vector<size_t> cluster;
26 
27 struct cluster_group {
28  std::vector<cluster> clusters;
30 };
31 
33  std::vector<triplet> triplets;
34  std::vector<cluster_group> history;
35 };
36 
38  float scale;
39 
40 public:
41  ScaleTripletMetric(float s);
42  float operator()(triplet const &lhs, triplet const &rhs, pcl::PointCloud<pcl::PointXYZI>::ConstPtr cloud);
43 };
44 
45 typedef float (*ClusterMetric)(cluster const &lhs, cluster const &rhs, Eigen::MatrixXf const &d,
46  pcl::PointCloud<pcl::PointXYZI>::ConstPtr cloud);
47 
48 inline float singleLinkClusterMetric(cluster const &lhs, cluster const &rhs, Eigen::MatrixXf const &d,
49  pcl::PointCloud<pcl::PointXYZI>::ConstPtr cloud)
50 {
51  float result = std::numeric_limits<float>::infinity();
52 
53  for (std::vector<size_t>::const_iterator a = lhs.begin(); a < lhs.end(); a++) {
54  for (std::vector<size_t>::const_iterator b = rhs.begin(); b < rhs.end(); b++) {
55  float distance = d(*a, *b);
56 
57  if (distance < result)
58  result = distance;
59  }
60  }
61 
62  return result;
63 }
64 
65 inline float completeLinkClusterMetric(cluster const &lhs, cluster const &rhs, Eigen::MatrixXf const &d,
66  pcl::PointCloud<pcl::PointXYZI>::ConstPtr cloud)
67 {
68  float result = 0.0f;
69 
70  for (std::vector<size_t>::const_iterator a = lhs.begin(); a < lhs.end(); a++) {
71  for (std::vector<size_t>::const_iterator b = rhs.begin(); b < rhs.end(); b++) {
72  float distance = d(*a, *b);
73 
74  if (distance > result)
75  result = distance;
76  }
77  }
78 
79  return result;
80 }
81 
82 std::vector<triplet> generateTriplets(pcl::PointCloud<pcl::PointXYZI>::ConstPtr cloud, size_t nnKandidates = 12,
83  size_t nBest = 2, float maxError = 0.015f);
84 cluster_group cleanupClusterGroup(cluster_group const &clusterGroup, size_t minTriplets = 7);
85 Cluster
86 toCluster(std::vector<hc::triplet> const &triplets, hc::cluster_group const &clusterGroup, size_t pointIndexCount);
87 cluster_group compute_hc(pcl::PointCloud<pcl::PointXYZI>::ConstPtr cloud, std::vector<triplet> const &triplets,
88  ScaleTripletMetric triplet_metric, float cdist, int opt_verbose);
89 } // namespace hc
90 
91 #endif
hc::triplet::pointIndexC
size_t pointIndexC
Definition: hc.h:18
hc::completeLinkClusterMetric
float completeLinkClusterMetric(cluster const &lhs, cluster const &rhs, Eigen::MatrixXf const &d, pcl::PointCloud< pcl::PointXYZI >::ConstPtr cloud)
Definition: hc.h:65
hc::triplet::pointIndexA
size_t pointIndexA
Definition: hc.h:16
f
double(* f)(double t, const double *par)
Definition: lmcurve.cxx:21
hc::cluster_history::history
std::vector< cluster_group > history
Definition: hc.h:34
hc
Definition: AtTrackFinderHC.h:32
hc::singleLinkClusterMetric
float singleLinkClusterMetric(cluster const &lhs, cluster const &rhs, Eigen::MatrixXf const &d, pcl::PointCloud< pcl::PointXYZI >::ConstPtr cloud)
Definition: hc.h:48
hc::cluster_group::clusters
std::vector< cluster > clusters
Definition: hc.h:28
hc::ClusterMetric
float(* ClusterMetric)(cluster const &lhs, cluster const &rhs, Eigen::MatrixXf const &d, pcl::PointCloud< pcl::PointXYZI >::ConstPtr cloud)
Definition: hc.h:45
hc::generateTriplets
std::vector< triplet > generateTriplets(pcl::PointCloud< pcl::PointXYZI >::ConstPtr cloud, size_t k, size_t n, float a)
Definition: hc.cxx:22
hc::triplet::direction
Eigen::Vector3f direction
Definition: hc.h:20
hc::cluster
std::vector< size_t > cluster
Definition: hc.h:25
Cluster
Definition: cluster.h:9
hc::triplet::error
float error
Definition: hc.h:21
hc::cleanupClusterGroup
cluster_group cleanupClusterGroup(cluster_group const &clusterGroup, size_t m)
Definition: hc.cxx:203
hc::cluster_history
Definition: hc.h:32
hc::ScaleTripletMetric::ScaleTripletMetric
ScaleTripletMetric(float s)
Definition: hc.cxx:245
hc::ScaleTripletMetric::operator()
float operator()(triplet const &lhs, triplet const &rhs, pcl::PointCloud< pcl::PointXYZI >::ConstPtr cloud)
Definition: hc.cxx:250
hc::triplet::center
Eigen::Vector3f center
Definition: hc.h:19
ScaleTripletMetric
Definition: triplet.h:32
hc::compute_hc
cluster_group compute_hc(pcl::PointCloud< pcl::PointXYZI >::ConstPtr cloud, std::vector< triplet > const &triplets, ScaleTripletMetric triplet_metric, float t, int opt_verbose)
Definition: hc.cxx:139
hc::cluster_group::bestClusterDistance
float bestClusterDistance
Definition: hc.h:29
hc::triplet::pointIndexB
size_t pointIndexB
Definition: hc.h:17
hc::cluster_history::triplets
std::vector< triplet > triplets
Definition: hc.h:33
hc::triplet::operator<
friend bool operator<(const triplet &t1, const triplet &t2)
Definition: hc.h:22
hc::toCluster
Cluster toCluster(std::vector< hc::triplet > const &triplets, hc::cluster_group const &clusterGroup, size_t pointIndexCount)
Definition: hc.cxx:216
hc::triplet
Definition: hc.h:15
cluster_group
std::vector< cluster_t > cluster_group
Definition: cluster.h:25
cluster.h
hc::cluster_group
Definition: hc.h:27
hc::ScaleTripletMetric
Definition: hc.h:37