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