ATTPCROOT  0.3.0-alpha
A ROOT-based framework for analyzing data from active target detectors
mst.cpp
Go to the documentation of this file.
1 #include "mst.h"
2 
3 #include <numeric>
4 
5 namespace mst
6 {
7  std::vector<edge> calculateSquareDistances(pcl::PointCloud<pcl::PointXYZI>::ConstPtr cloud, MstMetric metric)
8  {
9  const size_t cloud_size = cloud->size();
10  std::vector<edge> result;
11  result.reserve(cloud_size * (cloud_size - 1) / 2);
12 
13  for(size_t i = 0; i < cloud_size; ++i)
14  {
15  pcl::PointXYZI pointA = (*cloud)[i];
16 
17  for(size_t j = i + 1; j < cloud_size; ++j)
18  {
19  pcl::PointXYZI pointB = (*cloud)[j];
20 
21  float distance = metric(pointA, pointB, i, j, cloud);
22 
23  edge e = {
24  i, // voxelIndexA
25  j, // voxelIndexB
26  distance // distance
27  };
28 
29  result.push_back(e);
30  }
31  }
32 
33  std::sort(result.begin(), result.end());
34 
35  return result;
36  }
37 
38  std::vector<state> calculateMinimumSpanningTree(pcl::PointCloud<pcl::PointXYZI>::ConstPtr cloud, MstMetric metric, MstFilter filter)
39  {
40  std::vector<edge> edges = calculateSquareDistances(cloud, metric);
41  std::vector<state> result;
42 
43  std::vector<edge> selectedEdges;
44  std::vector<size_t> groups(cloud->size());
45  for(size_t i = 0; i < groups.size(); i++) {
46  groups[i] = i;
47  }
48 
49  for(std::vector<edge>::const_iterator edge = edges.begin(); edge != edges.end(); ++edge)
50  {
51  size_t groupA = groups[edge->voxelIndexA];
52  size_t groupB = groups[edge->voxelIndexB];
53 
54  // if no circle
55  if(groupA != groupB && filter(cloud, *edge, selectedEdges, groups))
56  {
57  selectedEdges.push_back(*edge);
58 
59  // merge groups
60  for(std::vector<size_t>::iterator it = groups.begin(); it != groups.end(); ++it)
61  {
62  if(*it == groupB)
63  *it = groupA;
64  }
65 
66  state state = {
67  selectedEdges, // edges
68  groups, // groups
69  };
70 
71  result.push_back(state);
72  }
73  }
74 
75  return result;
76  }
77 }
mst::state
Definition: mst.h:14
mst::MstFilter
bool(* MstFilter)(pcl::PointCloud< pcl::PointXYZI >::ConstPtr, edge const &, std::vector< edge > const &, std::vector< size_t > const &)
Definition: mst.h:21
mst::MstMetric
float(* MstMetric)(pcl::PointXYZI const &, pcl::PointXYZI const &, size_t, size_t, pcl::PointCloud< pcl::PointXYZI >::ConstPtr)
Definition: mst.h:19
mst::edge
Definition: mst.h:8
mst::edge::voxelIndexA
size_t voxelIndexA
Definition: mst.h:9
mst::calculateMinimumSpanningTree
std::vector< state > calculateMinimumSpanningTree(pcl::PointCloud< pcl::PointXYZI >::ConstPtr cloud, MstMetric metric, MstFilter filter)
Definition: mst.cpp:38
mst
Definition: mst.cpp:6
mst.h
mst::edge::voxelIndexB
size_t voxelIndexB
Definition: mst.h:9
mst::calculateSquareDistances
std::vector< edge > calculateSquareDistances(pcl::PointCloud< pcl::PointXYZI >::ConstPtr cloud, MstMetric metric)
Definition: mst.cpp:7