All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
ClusterPathFinder_tool.cc
Go to the documentation of this file.
1 /**
2  * @file ClusterPathFinder_tool.cc
3  *
4  * @brief art Tool for comparing clusters and merging those that are consistent
5  *
6  */
7 
8 // Framework Includes
9 #include "art/Utilities/ToolMacros.h"
10 #include "art/Utilities/make_tool.h"
11 #include "cetlib/cpu_timer.h"
12 #include "fhiclcpp/ParameterSet.h"
13 #include "messagefacility/MessageLogger/MessageLogger.h"
14 
18 
19 // LArSoft includes
22 
23 // Eigen
24 #include <Eigen/Core>
25 
26 // std includes
27 #include <string>
28 #include <iostream>
29 #include <memory>
30 
31 //------------------------------------------------------------------------------------------------------------------------------------------
32 // implementation follows
33 
34 namespace lar_cluster3d {
35 
36 class ClusterPathFinder : virtual public IClusterModAlg
37 {
38 public:
39  /**
40  * @brief Constructor
41  *
42  * @param pset
43  */
44  explicit ClusterPathFinder(const fhicl::ParameterSet&);
45 
46  /**
47  * @brief Destructor
48  */
50 
51  void configure(fhicl::ParameterSet const &pset) override;
52 
53  /**
54  * @brief Interface for initializing histograms if they are desired
55  * Note that the idea is to put hisgtograms in a subfolder
56  *
57  * @param TFileDirectory - the folder to store the hists in
58  */
59  void initializeHistograms(art::TFileDirectory&) override;
60 
61  /**
62  * @brief Scan an input collection of clusters and modify those according
63  * to the specific implementing algorithm
64  *
65  * @param clusterParametersList A list of cluster objects (parameters from associated hits)
66  */
67  void ModifyClusters(reco::ClusterParametersList&) const override;
68 
69  /**
70  * @brief If monitoring, recover the time to execute a particular function
71  */
72  float getTimeToExecute() const override {return m_timeToProcess;}
73 
74 private:
75 
76  /**
77  * @brief Use PCA to try to find path in cluster
78  *
79  * @param clusterParameters The given cluster parameters object to try to split
80  * @param clusterParametersList The list of clusters
81  */
82  reco::ClusterParametersList::iterator breakIntoTinyBits(reco::ClusterParameters& cluster,
83  reco::ClusterParametersList::iterator positionItr,
84  reco::ClusterParametersList& outputClusterList,
85  int level = 0) const;
86 
87  float closestApproach(const Eigen::Vector3f&,
88  const Eigen::Vector3f&,
89  const Eigen::Vector3f&,
90  const Eigen::Vector3f&,
91  Eigen::Vector3f&,
92  Eigen::Vector3f&) const;
93 
94  using Point = std::tuple<float,float,const reco::ClusterHit3D*>;
95  using PointList = std::list<Point>;
96  using MinMaxPoints = std::pair<Point,Point>;
97  using MinMaxPointPair = std::pair<MinMaxPoints,MinMaxPoints>;
98 
99  void buildConvexHull(reco::ClusterParameters& clusterParameters, int level = 0) const;
100  void buildVoronoiDiagram(reco::ClusterParameters& clusterParameters, int level = 0) const;
101  /**
102  * @brief Data members to follow
103  */
105  size_t m_minTinyClusterSize; ///< Minimum size for a "tiny" cluster
106  mutable float m_timeToProcess; ///<
107 
108  std::unique_ptr<lar_cluster3d::IClusterAlg> m_clusterAlg; ///< Algorithm to do 3D space point clustering
109  PrincipalComponentsAlg m_pcaAlg; // For running Principal Components Analysis
110 };
111 
112 ClusterPathFinder::ClusterPathFinder(fhicl::ParameterSet const &pset) :
113  m_pcaAlg(pset.get<fhicl::ParameterSet>("PrincipalComponentsAlg"))
114 {
115  this->configure(pset);
116 }
117 
118 //------------------------------------------------------------------------------------------------------------------------------------------
119 
121 {
122 }
123 
124 //------------------------------------------------------------------------------------------------------------------------------------------
125 
126 void ClusterPathFinder::configure(fhicl::ParameterSet const &pset)
127 {
128  m_enableMonitoring = pset.get<bool> ("EnableMonitoring", true );
129  m_minTinyClusterSize = pset.get<size_t>("MinTinyClusterSize",40);
130  m_clusterAlg = art::make_tool<lar_cluster3d::IClusterAlg>(pset.get<fhicl::ParameterSet>("ClusterAlg"));
131 
132  m_timeToProcess = 0.;
133 
134  return;
135 }
136 
137 void ClusterPathFinder::initializeHistograms(art::TFileDirectory&)
138 {
139  return;
140 }
141 
143 {
144  /**
145  * @brief Top level interface for algorithm to consider pairs of clusters from the input
146  * list and determine if they are consistent with each other and, therefore, should
147  * be merged. This is done by looking at the PCA for each cluster and looking at the
148  * projection of the primary axis along the vector connecting their centers.
149  */
150 
151  // Initial clustering is done, now trim the list and get output parameters
152  cet::cpu_timer theClockBuildClusters;
153 
154  // Start clocks if requested
155  if (m_enableMonitoring) theClockBuildClusters.start();
156 
157  int countClusters(0);
158 
159  // This is the loop over candidate 3D clusters
160  // Note that it might be that the list of candidate clusters is modified by splitting
161  // So we use the following construct to make sure we get all of them
162  reco::ClusterParametersList::iterator clusterParametersListItr = clusterParametersList.begin();
163 
164  while(clusterParametersListItr != clusterParametersList.end())
165  {
166  // Dereference to get the cluster paramters
167  reco::ClusterParameters& clusterParameters = *clusterParametersListItr;
168 
169  std::cout << "**> Looking at Cluster " << countClusters++ << std::endl;
170 
171  // It turns out that computing the convex hull surrounding the points in the 2D projection onto the
172  // plane of largest spread in the PCA is a good way to break up the cluster... and we do it here since
173  // we (currently) want this to be part of the standard output
174  buildVoronoiDiagram(clusterParameters);
175 
176  // Make sure our cluster has enough hits...
177  if (clusterParameters.getHitPairListPtr().size() > m_minTinyClusterSize)
178  {
179  // Get an interim cluster list
180  reco::ClusterParametersList reclusteredParameters;
181 
182  // Call the main workhorse algorithm for building the local version of candidate 3D clusters
183  m_clusterAlg->Cluster3DHits(clusterParameters.getHitPairListPtr(), reclusteredParameters);
184 
185  std::cout << ">>>>>>>>>>> Reclustered to " << reclusteredParameters.size() << " Clusters <<<<<<<<<<<<<<<" << std::endl;
186 
187  // Only process non-empty results
188  if (!reclusteredParameters.empty())
189  {
190  // Loop over the reclustered set
191  for (auto& cluster : reclusteredParameters)
192  {
193  std::cout << "****> Calling breakIntoTinyBits" << std::endl;
194 
195  // Break our cluster into smaller elements...
196  breakIntoTinyBits(cluster, cluster.daughterList().end(), cluster.daughterList(), 4);
197 
198  std::cout << "****> Broke Cluster with " << cluster.getHitPairListPtr().size() << " into " << cluster.daughterList().size() << " sub clusters";
199  for(auto& clus : cluster.daughterList()) std::cout << ", " << clus.getHitPairListPtr().size();
200  std::cout << std::endl;
201 
202  // Add the daughters to the cluster
203  clusterParameters.daughterList().insert(clusterParameters.daughterList().end(),cluster);
204  }
205  }
206  }
207 
208  // Go to next cluster parameters object
209  clusterParametersListItr++;
210  }
211 
212  if (m_enableMonitoring)
213  {
214  theClockBuildClusters.stop();
215 
216  m_timeToProcess = theClockBuildClusters.accumulated_real_time();
217  }
218 
219  mf::LogDebug("Cluster3D") << ">>>>> Cluster Path finding done" << std::endl;
220 
221  return;
222 }
223 
224 reco::ClusterParametersList::iterator ClusterPathFinder::breakIntoTinyBits(reco::ClusterParameters& clusterToBreak,
225  reco::ClusterParametersList::iterator positionItr,
226  reco::ClusterParametersList& outputClusterList,
227  int level) const
228 {
229  // This needs to be a recursive routine...
230  // Idea is to take the input cluster and order 3D hits by arclength along PCA primary axis
231  // If the cluster is above the minimum number of hits then we divide into two and call ourself
232  // with the two halves. This means we form a new cluster with hits and PCA and then call ourself
233  // If the cluster is below the minimum then we can't break any more, simply add this cluster to
234  // the new output list.
235 
236  // set an indention string
237  std::string pluses(level/2, '+');
238  std::string indent(level/2, ' ');
239 
240  indent += pluses;
241 
242  reco::ClusterParametersList::iterator inputPositionItr = positionItr;
243 
244  // To make best use of this we'll also want the PCA for this cluster... so...
245  // Recover the prime ingredients
246  reco::PrincipalComponents& fullPCA = clusterToBreak.getFullPCA();
247 
248  std::cout << indent << ">>> breakIntoTinyBits with " << clusterToBreak.getHitPairListPtr().size() << " input hits " << std::endl;
249 
250  // It turns out that computing the convex hull surrounding the points in the 2D projection onto the
251  // plane of largest spread in the PCA is a good way to break up the cluster... and we do it here since
252  // we (currently) want this to be part of the standard output
253  buildConvexHull(clusterToBreak, level+2);
254 
255  bool storeCurrentCluster(true);
256  int minimumClusterSize(m_minTinyClusterSize);
257 
258  // Create a rough cut intended to tell us when we've reached the land of diminishing returns
259  if (clusterToBreak.getBestEdgeList().size() > 6 &&
260  clusterToBreak.getHitPairListPtr().size() > size_t(2 * minimumClusterSize))
261  {
262  // Recover the list of 3D hits associated to this cluster
263  reco::HitPairListPtr& clusHitPairVector = clusterToBreak.getHitPairListPtr();
264 
265  // Calculate the doca to the PCA primary axis for each 3D hit
266  // Importantly, this also gives us the arclength along that axis to the hit
267  m_pcaAlg.PCAAnalysis_calc3DDocas(clusHitPairVector, fullPCA);
268 
269  // Sort the hits along the PCA
270  clusHitPairVector.sort([](const auto& left, const auto& right){return left->getArclenToPoca() < right->getArclenToPoca();});
271 
272  // Now we use the convex hull vertex points to form split points for breaking up the incoming cluster
273  reco::EdgeList& bestEdgeList = clusterToBreak.getBestEdgeList();
274  std::vector<const reco::ClusterHit3D*> vertexHitVec;
275 
276  std::cout << indent << "+> Breaking cluster, convex hull has " << bestEdgeList.size() << " edges to work with" << std::endl;
277 
278  for(const auto& edge : bestEdgeList)
279  {
280  vertexHitVec.push_back(std::get<0>(edge));
281  vertexHitVec.push_back(std::get<1>(edge));
282  }
283 
284  // Sort this vector, we aren't worried about duplicates right now...
285  std::sort(vertexHitVec.begin(),vertexHitVec.end(),[](const auto& left, const auto& right){return left->getArclenToPoca() < right->getArclenToPoca();});
286 
287  // Now we create a list of pairs of iterators to the start and end of each subcluster
288  using Hit3DItrPair = std::pair<reco::HitPairListPtr::iterator,reco::HitPairListPtr::iterator>;
289  using VertexPairList = std::list<Hit3DItrPair>;
290 
291  VertexPairList vertexPairList;
292  reco::HitPairListPtr::iterator firstHitItr = clusHitPairVector.begin();
293 
294  for(const auto& hit3D : vertexHitVec)
295  {
296  reco::HitPairListPtr::iterator vertexItr = std::find(firstHitItr,clusHitPairVector.end(),hit3D);
297 
298  if (vertexItr == clusHitPairVector.end())
299  {
300  std::cout << indent << ">>>>>>>>>>>>>>>>> Hit not found in input list, cannot happen? <<<<<<<<<<<<<<<<<<<" << std::endl;
301  break;
302  }
303 
304  std::cout << indent << "+> -- Distance from first to current vertex point: " << std::distance(firstHitItr,vertexItr) << " first: " << *firstHitItr << ", vertex: " << *vertexItr;
305 
306  // Require a minimum number of points...
307  if (std::distance(firstHitItr,vertexItr) > minimumClusterSize)
308  {
309  vertexPairList.emplace_back(Hit3DItrPair(firstHitItr,vertexItr));
310  firstHitItr = vertexItr;
311 
312  std::cout << " ++ made pair ";
313  }
314 
315  std::cout << std::endl;
316  }
317 
318  // Not done if there is distance from first to end of list
319  if (std::distance(firstHitItr,clusHitPairVector.end()) > 0)
320  {
321  std::cout << indent << "+> loop over vertices done, remant distance: " << std::distance(firstHitItr,clusHitPairVector.end()) << std::endl;
322 
323  // In the event we don't have the minimum number of hits we simply extend the last pair
324  if (!vertexPairList.empty() && std::distance(firstHitItr,clusHitPairVector.end()) < minimumClusterSize)
325  vertexPairList.back().second = clusHitPairVector.end();
326  else
327  vertexPairList.emplace_back(Hit3DItrPair(firstHitItr,clusHitPairVector.end()));
328  }
329 
330  std::cout << indent << "+> ---> breaking cluster into " << vertexPairList.size() << " subclusters" << std::endl;
331 
332  if (vertexPairList.size() > 1)
333  {
334  storeCurrentCluster = false;
335 
336  // Ok, now loop through our pairs
337  for(auto& hit3DItrPair : vertexPairList)
338  {
339  reco::ClusterParameters clusterParams;
340  reco::HitPairListPtr& hitPairListPtr = clusterParams.getHitPairListPtr();
341 
342  std::cout << indent << "+> -- building new cluster, size: " << std::distance(hit3DItrPair.first,hit3DItrPair.second) << std::endl;
343 
344  // size the container...
345  hitPairListPtr.resize(std::distance(hit3DItrPair.first,hit3DItrPair.second));
346 
347  // and copy the hits into the container
348  std::copy(hit3DItrPair.first,hit3DItrPair.second,hitPairListPtr.begin());
349 
350  // First stage of feature extraction runs here
351  m_pcaAlg.PCAAnalysis_3D(hitPairListPtr, clusterParams.getFullPCA());
352 
353  // Recover the new fullPCA
354  reco::PrincipalComponents& newFullPCA = clusterParams.getFullPCA();
355 
356  // Must have a valid pca
357  if (newFullPCA.getSvdOK())
358  {
359  std::cout << indent << "+> -- >> cluster has a valid Full PCA" << std::endl;
360 
361  // Need to check if the PCA direction has been reversed
362  Eigen::Vector3f fullPrimaryVec(fullPCA.getEigenVectors().row(0));
363  Eigen::Vector3f newPrimaryVec(newFullPCA.getEigenVectors().row(0));
364 
365  // If the PCA's are opposite the flip the axes
366  if (fullPrimaryVec.dot(newPrimaryVec) < 0.)
367  {
368  for(size_t vecIdx = 0; vecIdx < 3; vecIdx++) newFullPCA.flipAxis(vecIdx);
369  }
370 
371  // Set the skeleton PCA to make sure it has some value
372  clusterParams.getSkeletonPCA() = clusterParams.getFullPCA();
373 
374  positionItr = breakIntoTinyBits(clusterParams, positionItr, outputClusterList, level+4);
375  }
376  }
377  }
378  }
379 
380  // First question, are we done breaking?
381  if (storeCurrentCluster)
382  {
383  // I think this is where we fill out the rest of the parameters?
384  // Start by adding the 2D hits...
385  // See if we can avoid duplicates by temporarily transferring to a set
386  std::set<const reco::ClusterHit2D*> hitSet;
387 
388  // Loop through 3D hits to get a set of unique 2D hits
389  for(const auto& hit3D : clusterToBreak.getHitPairListPtr())
390  {
391  for(const auto& hit2D : hit3D->getHits())
392  {
393  if (hit2D) hitSet.insert(hit2D);
394  }
395  }
396 
397  // Now add these to the new cluster
398  for(const auto& hit2D : hitSet)
399  {
400  hit2D->setStatusBit(reco::ClusterHit2D::USED);
401  clusterToBreak.UpdateParameters(hit2D);
402  }
403 
404  std::cout << indent << "*********>>> storing new subcluster of size " << clusterToBreak.getHitPairListPtr().size() << std::endl;
405 
406  positionItr = outputClusterList.insert(positionItr,clusterToBreak);
407 
408  // The above points to the element, want the next element
409  positionItr++;
410  }
411  else if (inputPositionItr == positionItr) std::cout << indent << "***** DID NOT STORE A CLUSTER *****" << std::endl;
412 
413  return positionItr;
414 }
415 
416 void ClusterPathFinder::buildConvexHull(reco::ClusterParameters& clusterParameters, int level) const
417 {
418  // set an indention string
419  std::string minuses(level/2, '-');
420  std::string indent(level/2, ' ');
421 
422  indent += minuses;
423 
424  // The plan is to build the enclosing 2D polygon around the points in the PCA plane of most spread for this cluster
425  // To do so we need to start by building a list of 2D projections onto the plane of most spread...
426  reco::PrincipalComponents& pca = clusterParameters.getFullPCA();
427 
428  // Recover the parameters from the Principal Components Analysis that we need to project and accumulate
429  const Eigen::Vector3f& pcaCenter = pca.getAvePosition();
430 
431  //dcel2d::PointList pointList;
432  using Point = std::tuple<float,float,const reco::ClusterHit3D*>;
433  using PointList = std::list<Point>;
434  PointList pointList;
435 
436  // Loop through hits and do projection to plane
437  for(const auto& hit3D : clusterParameters.getHitPairListPtr())
438  {
439  Eigen::Vector3f pcaToHitVec(hit3D->getPosition()[0] - pcaCenter(0),
440  hit3D->getPosition()[1] - pcaCenter(1),
441  hit3D->getPosition()[2] - pcaCenter(2));
442  Eigen::Vector3f pcaToHit = pca.getEigenVectors() * pcaToHitVec;
443 
444  pointList.emplace_back(dcel2d::Point(pcaToHit(1),pcaToHit(2),hit3D));
445  }
446 
447  // Sort the point vec by increasing x, then increase y
448  pointList.sort([](const auto& left, const auto& right){return (std::abs(std::get<0>(left) - std::get<0>(right)) > std::numeric_limits<float>::epsilon()) ? std::get<0>(left) < std::get<0>(right) : std::get<1>(left) < std::get<1>(right);});
449 
450  // containers for finding the "best" hull...
451  std::vector<ConvexHull> convexHullVec;
452  std::vector<PointList> rejectedListVec;
453  bool increaseDepth(pointList.size() > 5);
454  float lastArea(std::numeric_limits<float>::max());
455 
456  while(increaseDepth)
457  {
458  // Get another convexHull container
459  convexHullVec.push_back(ConvexHull(pointList));
460  rejectedListVec.push_back(PointList());
461 
462  const ConvexHull& convexHull = convexHullVec.back();
463  PointList& rejectedList = rejectedListVec.back();
464  const PointList& convexHullPoints = convexHull.getConvexHull();
465 
466  increaseDepth = false;
467 
468  if (convexHull.getConvexHullArea() > 0.)
469  {
470  std::cout << indent << "-> built convex hull, 3D hits: " << pointList.size() << " with " << convexHullPoints.size() << " vertices" << ", area: " << convexHull.getConvexHullArea() << std::endl;
471  std::cout << indent << "-> -Points:";
472  for(const auto& point : convexHullPoints)
473  std::cout << " (" << std::get<0>(point) << "," << std::get<1>(point) << ")";
474  std::cout << std::endl;
475 
476  if (convexHullVec.size() < 2 || convexHull.getConvexHullArea() < 0.8 * lastArea)
477  {
478  for(auto& point : convexHullPoints)
479  {
480  pointList.remove(point);
481  rejectedList.emplace_back(point);
482  }
483  lastArea = convexHull.getConvexHullArea();
484 // increaseDepth = true;
485  }
486  }
487  }
488 
489  // do we have a valid convexHull?
490  while(!convexHullVec.empty() && convexHullVec.back().getConvexHullArea() < 0.5)
491  {
492  convexHullVec.pop_back();
493  rejectedListVec.pop_back();
494  }
495 
496  // If we found the convex hull then build edges around the region
497  if (!convexHullVec.empty())
498  {
499  size_t nRejectedTotal(0);
500  reco::HitPairListPtr hitPairListPtr = clusterParameters.getHitPairListPtr();
501 
502  for(const auto& rejectedList : rejectedListVec)
503  {
504  nRejectedTotal += rejectedList.size();
505 
506  for(const auto& rejectedPoint : rejectedList)
507  {
508  std::cout << indent << "-> -- Point is " << convexHullVec.back().findNearestDistance(rejectedPoint) << " from nearest edge" << std::endl;
509 
510  if (convexHullVec.back().findNearestDistance(rejectedPoint) > 0.5)
511  hitPairListPtr.remove(std::get<2>(rejectedPoint));
512  }
513  }
514 
515  std::cout << indent << "-> Removed " << nRejectedTotal << " leaving " << pointList.size() << "/" << hitPairListPtr.size() << " points" << std::endl;
516 
517  // Now add "edges" to the cluster to describe the convex hull (for the display)
518  reco::Hit3DToEdgeMap& edgeMap = clusterParameters.getHit3DToEdgeMap();
519  reco::EdgeList& bestEdgeList = clusterParameters.getBestEdgeList();
520 
521  Point lastPoint = convexHullVec.back().getConvexHull().front();
522 
523  for(auto& curPoint : convexHullVec.back().getConvexHull())
524  {
525  if (curPoint == lastPoint) continue;
526 
527  const reco::ClusterHit3D* lastPoint3D = std::get<2>(lastPoint);
528  const reco::ClusterHit3D* curPoint3D = std::get<2>(curPoint);
529 
530  float distBetweenPoints = (curPoint3D->getPosition()[0] - lastPoint3D->getPosition()[0]) * (curPoint3D->getPosition()[0] - lastPoint3D->getPosition()[0])
531  + (curPoint3D->getPosition()[1] - lastPoint3D->getPosition()[1]) * (curPoint3D->getPosition()[1] - lastPoint3D->getPosition()[1])
532  + (curPoint3D->getPosition()[2] - lastPoint3D->getPosition()[2]) * (curPoint3D->getPosition()[2] - lastPoint3D->getPosition()[2]);
533 
534  distBetweenPoints = std::sqrt(distBetweenPoints);
535 
536  reco::EdgeTuple edge(lastPoint3D,curPoint3D,distBetweenPoints);
537 
538  edgeMap[lastPoint3D].push_back(edge);
539  edgeMap[curPoint3D].push_back(edge);
540  bestEdgeList.emplace_back(edge);
541 
542  lastPoint = curPoint;
543  }
544  }
545 
546  return;
547 }
548 
549 
550 void ClusterPathFinder::buildVoronoiDiagram(reco::ClusterParameters& clusterParameters, int level) const
551 {
552  // The plan is to build the enclosing 2D polygon around the points in the PCA plane of most spread for this cluster
553  // To do so we need to start by building a list of 2D projections onto the plane of most spread...
554  reco::PrincipalComponents& pca = clusterParameters.getFullPCA();
555 
556  // Recover the parameters from the Principal Components Analysis that we need to project and accumulate
557  const Eigen::Vector3f& pcaCenter = pca.getAvePosition();
558  dcel2d::PointList pointList;
559 
560  // Loop through hits and do projection to plane
561  for(const auto& hit3D : clusterParameters.getHitPairListPtr())
562  {
563  Eigen::Vector3f pcaToHitVec(hit3D->getPosition()[0] - pcaCenter(0),
564  hit3D->getPosition()[1] - pcaCenter(1),
565  hit3D->getPosition()[2] - pcaCenter(2));
566  Eigen::Vector3f pcaToHit = pca.getEigenVectors() * pcaToHitVec;
567 
568  pointList.emplace_back(dcel2d::Point(pcaToHit(1),pcaToHit(2),hit3D));
569  }
570 
571  // Sort the point vec by increasing x, then increase y
572  pointList.sort([](const auto& left, const auto& right){return (std::abs(std::get<0>(left) - std::get<0>(right)) > std::numeric_limits<float>::epsilon()) ? std::get<0>(left) < std::get<0>(right) : std::get<1>(left) < std::get<1>(right);});
573 
574  // Set up the voronoi diagram builder
575  voronoi2d::VoronoiDiagram voronoiDiagram(clusterParameters.getHalfEdgeList(),clusterParameters.getVertexList(),clusterParameters.getFaceList());
576 
577  // And make the diagram
578  voronoiDiagram.buildVoronoiDiagram(pointList);
579 
580  // Recover the voronoi diagram vertex list and the container to store them in
581  dcel2d::VertexList& vertexList = clusterParameters.getVertexList();
582 
583  // Now get the inverse of the rotation matrix so we can get the vertex positions,
584  // which lie in the plane of the two largest PCA axes, in the standard coordinate system
585  Eigen::Matrix3f rotationMatrixInv = pca.getEigenVectors().inverse();
586 
587  // Translate and fill
588  for(auto& vertex : vertexList)
589  {
590  Eigen::Vector3f coords = rotationMatrixInv * vertex.getCoords();
591 
592  coords += pcaCenter;
593 
594  vertex.setCoords(coords);
595  }
596 
597  // Now do the Convex Hull
598  // Now add "edges" to the cluster to describe the convex hull (for the display)
599  reco::Hit3DToEdgeMap& edgeMap = clusterParameters.getHit3DToEdgeMap();
600  reco::EdgeList& bestEdgeList = clusterParameters.getBestEdgeList();
601 
602 // const dcel2d::PointList& edgePoints = voronoiDiagram.getConvexHull();
603 // PointList localList;
604 //
605 // for(const auto& edgePoint : edgePoints) localList.emplace_back(std::get<0>(edgePoint),std::get<1>(edgePoint),std::get<2>(edgePoint));
606 //
607 // // Sort the point vec by increasing x, then increase y
608 // localList.sort([](const auto& left, const auto& right){return (std::abs(std::get<0>(left) - std::get<0>(right)) > std::numeric_limits<float>::epsilon()) ? std::get<0>(left) < std::get<0>(right) : std::get<1>(left) < std::get<1>(right);});
609 //
610 // // Why do I need to do this?
611 // ConvexHull convexHull(localList);
612 //
613 // Point lastPoint = convexHull.getConvexHull().front();
614  dcel2d::Point lastPoint = voronoiDiagram.getConvexHull().front();
615 
616 // std::cout << "@@@@>> Build convex hull, voronoi handed " << edgePoints.size() << " points, convexHull cut to " << convexHull.getConvexHull().size() << std::endl;
617 
618 // for(auto& curPoint : convexHull.getConvexHull())
619  for(auto& curPoint : voronoiDiagram.getConvexHull())
620  {
621  if (curPoint == lastPoint) continue;
622 
623  const reco::ClusterHit3D* lastPoint3D = std::get<2>(lastPoint);
624  const reco::ClusterHit3D* curPoint3D = std::get<2>(curPoint);
625 
626  float distBetweenPoints = (curPoint3D->getPosition()[0] - lastPoint3D->getPosition()[0]) * (curPoint3D->getPosition()[0] - lastPoint3D->getPosition()[0])
627  + (curPoint3D->getPosition()[1] - lastPoint3D->getPosition()[1]) * (curPoint3D->getPosition()[1] - lastPoint3D->getPosition()[1])
628  + (curPoint3D->getPosition()[2] - lastPoint3D->getPosition()[2]) * (curPoint3D->getPosition()[2] - lastPoint3D->getPosition()[2]);
629 
630  distBetweenPoints = std::sqrt(distBetweenPoints);
631 
632  reco::EdgeTuple edge(lastPoint3D,curPoint3D,distBetweenPoints);
633 
634  edgeMap[lastPoint3D].push_back(edge);
635  edgeMap[curPoint3D].push_back(edge);
636  bestEdgeList.emplace_back(edge);
637 
638  lastPoint = curPoint;
639  }
640 
641  std::cout << "****> vertexList containted " << vertexList.size() << " vertices" << std::endl;
642 
643  return;
644 }
645 
646 
647 float ClusterPathFinder::closestApproach(const Eigen::Vector3f& P0,
648  const Eigen::Vector3f& u0,
649  const Eigen::Vector3f& P1,
650  const Eigen::Vector3f& u1,
651  Eigen::Vector3f& poca0,
652  Eigen::Vector3f& poca1) const
653 {
654  // Technique is to compute the arclength to each point of closest approach
655  Eigen::Vector3f w0 = P0 - P1;
656  float a(1.);
657  float b(u0.dot(u1));
658  float c(1.);
659  float d(u0.dot(w0));
660  float e(u1.dot(w0));
661  float den(a * c - b * b);
662 
663  float arcLen0 = (b * e - c * d) / den;
664  float arcLen1 = (a * e - b * d) / den;
665 
666  poca0 = P0 + arcLen0 * u0;
667  poca1 = P1 + arcLen1 * u1;
668 
669  return (poca0 - poca1).norm();
670 }
671 
672 DEFINE_ART_CLASS_TOOL(ClusterPathFinder)
673 } // namespace lar_cluster3d
process_name vertex
Definition: cheaterreco.fcl:51
void PCAAnalysis_calc3DDocas(const reco::HitPairListPtr &hitPairVector, const reco::PrincipalComponents &pca) const
bool getSvdOK() const
Definition: Cluster3D.h:243
void PCAAnalysis_3D(const reco::HitPairListPtr &hitPairList, reco::PrincipalComponents &pca, bool skeletonOnly=false) const
process_name cluster
Definition: cheaterreco.fcl:51
void flipAxis(size_t axis)
Definition: Cluster3D.cxx:208
const Eigen::Vector3f getPosition() const
Definition: Cluster3D.h:157
walls no right
Definition: selectors.fcl:105
void configure(fhicl::ParameterSet const &pset) override
Interface for configuring the particular algorithm tool.
reco::PrincipalComponents & getSkeletonPCA()
Definition: Cluster3D.h:477
std::unique_ptr< lar_cluster3d::IClusterAlg > m_clusterAlg
Algorithm to do 3D space point clustering.
reco::EdgeList & getBestEdgeList()
Definition: Cluster3D.h:480
float closestApproach(const Eigen::Vector3f &, const Eigen::Vector3f &, const Eigen::Vector3f &, const Eigen::Vector3f &, Eigen::Vector3f &, Eigen::Vector3f &) const
reco::Hit3DToEdgeMap & getHit3DToEdgeMap()
Definition: Cluster3D.h:478
reco::HitPairListPtr & getHitPairListPtr()
Definition: Cluster3D.h:475
std::tuple< float, float, const reco::ClusterHit3D * > Point
std::pair< MinMaxPoints, MinMaxPoints > MinMaxPointPair
IClusterModAlg interface class definiton.
process_name gaushit a
ClusterParametersList & daughterList()
Definition: Cluster3D.h:449
Implements a ConvexHull for use in clustering.
T abs(T value)
std::list< EdgeTuple > EdgeList
Definition: Cluster3D.h:344
void buildConvexHull(reco::ClusterParameters &clusterParameters, int level=0) const
void initializeHistograms(art::TFileDirectory &) override
Interface for initializing histograms if they are desired Note that the idea is to put hisgtograms in...
VoronoiDiagram class definiton.
Definition: Voronoi.h:31
double distance(geo::Point_t const &point, CathodeDesc_t const &cathode)
Returns the distance of a point from the cathode.
reco::PrincipalComponents & getFullPCA()
Definition: Cluster3D.h:476
float getTimeToExecute() const override
If monitoring, recover the time to execute a particular function.
std::tuple< double, double, const reco::ClusterHit3D * > Point
Definitions used by the VoronoiDiagram algorithm.
Definition: DCEL.h:44
void buildVoronoiDiagram(const dcel2d::PointList &)
Given an input set of 2D points construct a 2D voronoi diagram.
Definition: Voronoi.cxx:133
const Eigen::Vector3f & getAvePosition() const
Definition: Cluster3D.h:247
size_t m_minTinyClusterSize
Minimum size for a &quot;tiny&quot; cluster.
std::tuple< const reco::ClusterHit3D *, const reco::ClusterHit3D *, double > EdgeTuple
Definition: Cluster3D.h:343
std::list< const reco::ClusterHit3D * > HitPairListPtr
Definition: Cluster3D.h:334
void UpdateParameters(const reco::ClusterHit2D *hit)
Definition: Cluster3D.h:451
walls no left
Definition: selectors.fcl:105
float getConvexHullArea() const
recover the area of the convex hull
Definition: ConvexHull.h:78
This provides an art tool interface definition for 3D Cluster algorithms.
const PointList & getConvexHull() const
recover the list of convex hull vertices
Definition: ConvexHull.h:58
reco::ClusterParametersList::iterator breakIntoTinyBits(reco::ClusterParameters &cluster, reco::ClusterParametersList::iterator positionItr, reco::ClusterParametersList &outputClusterList, int level=0) const
Use PCA to try to find path in cluster.
std::pair< Point, Point > MinMaxPoints
ConvexHull class definiton.
Definition: ConvexHull.h:27
This header file defines the interface to a principal components analysis designed to be used within ...
bool m_enableMonitoring
Data members to follow.
void buildVoronoiDiagram(reco::ClusterParameters &clusterParameters, int level=0) const
dcel2d::FaceList & getFaceList()
Definition: Cluster3D.h:482
This provides an art tool interface definition for 3D Cluster algorithms.
do i e
T copy(T const &v)
void ModifyClusters(reco::ClusterParametersList &) const override
Scan an input collection of clusters and modify those according to the specific implementing algorith...
std::list< Point > PointList
Definition: DCEL.h:45
std::unordered_map< const reco::ClusterHit3D *, reco::EdgeList > Hit3DToEdgeMap
Definition: Cluster3D.h:346
process_name physics producers generator hPHist_pi physics producers generator P0
ClusterPathFinder(const fhicl::ParameterSet &)
Constructor.
dcel2d::HalfEdgeList & getHalfEdgeList()
Definition: Cluster3D.h:484
std::list< Vertex > VertexList
Definition: DCEL.h:182
BEGIN_PROLOG don t mess with this pandoraTrackGausCryoW true
BEGIN_PROLOG could also be cout
std::list< ClusterParameters > ClusterParametersList
Definition: Cluster3D.h:403
const EigenVectors & getEigenVectors() const
Definition: Cluster3D.h:246
dcel2d::VertexList & getVertexList()
Definition: Cluster3D.h:483