LArSoft  v06_85_00
Liquid Argon Software toolkit - http://larsoft.org/
VoronoiPathFinder_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 VoronoiPathFinder : virtual public IClusterModAlg
45 {
46 public:
52  explicit VoronoiPathFinder(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  // The plan we think we need here:
267  // 1) recover the list of edges
268  // 2) sort by length
269  // 3) then find the defect point with the furthest distance perpendicular to this edge
270  // 4) Use this as the breakpoint, sort points according to their position along the PCA and
271  // side of the perpendicular edge...
272  reco::EdgeList& bestEdgeList = clusterToBreak.getBestEdgeList();
273 
274  const reco::EdgeTuple& longEdge = *std::max_element(bestEdgeList.begin(),bestEdgeList.end(),[](const auto& first, const auto& second){return std::get<2>(first) < std::get<2>(second);});
275 
276  // Get a vector representing this edge
277  const reco::ClusterHit3D* firstEdgeHit = std::get<0>(longEdge);
278  const reco::ClusterHit3D* secondEdgeHit = std::get<1>(longEdge);
279  double edgeLen = std::get<2>(longEdge);
280  Eigen::Vector3f edgeVec(secondEdgeHit->getPosition()[0] - firstEdgeHit->getPosition()[0],
281  secondEdgeHit->getPosition()[1] - firstEdgeHit->getPosition()[1],
282  secondEdgeHit->getPosition()[2] - firstEdgeHit->getPosition()[2]);
283 
284  // normalize it
285  edgeVec.normalize();
286 
287  // Keep track of the winner
288  const reco::ClusterHit3D* furthestHit(NULL);
289  float furthestDistance(0.);
290 
291  // Now loop through all the edges and search for the furthers point
292  for(const auto& edge : bestEdgeList)
293  {
294  const reco::ClusterHit3D* nextEdgeHit = std::get<0>(edge); // recover the first point
295 
296  // Create vector to this point from the longest edge
297  Eigen::Vector3f hitToEdgeVec(nextEdgeHit->getPosition()[0] - firstEdgeHit->getPosition()[0],
298  nextEdgeHit->getPosition()[1] - firstEdgeHit->getPosition()[1],
299  nextEdgeHit->getPosition()[2] - firstEdgeHit->getPosition()[2]);
300 
301  // Get projection
302  float hitProjection = hitToEdgeVec.dot(edgeVec);
303 
304  // Require that the point is really "opposite" the longest edge
305  if (hitProjection > 0. && hitProjection < edgeLen)
306  {
307  Eigen::Vector3f distToHitVec = hitToEdgeVec - hitProjection * edgeVec;
308  float distToHit = distToHitVec.norm();
309 
310  if (distToHit > furthestDistance)
311  {
312  furthestDistance = distToHit;
313  furthestHit = nextEdgeHit;
314  }
315  }
316  }
317 
318  // Make sure we have a hit
319  if (furthestHit)
320  {
321  // Recover the list of 3D hits associated to this cluster
322  reco::HitPairListPtr& clusHitPairVector = clusterToBreak.getHitPairListPtr();
323 
324  // Calculate the doca to the PCA primary axis for each 3D hit
325  // Importantly, this also gives us the arclength along that axis to the hit
326  m_pcaAlg.PCAAnalysis_calc3DDocas(clusHitPairVector, fullPCA);
327 
328  // Sort the hits along the PCA
329  clusHitPairVector.sort([](const auto& left, const auto& right){return left->getArclenToPoca() < right->getArclenToPoca();});
330 
331  // Now find the hit identified above as furthest away
332  reco::HitPairListPtr::iterator vertexItr = std::find(clusHitPairVector.begin(),clusHitPairVector.end(),furthestHit);
333 
334  // Now we create a list of pairs of iterators to the start and end of each subcluster
335  using Hit3DItrPair = std::pair<reco::HitPairListPtr::iterator,reco::HitPairListPtr::iterator>;
336  using VertexPairList = std::list<Hit3DItrPair>;
337 
338  VertexPairList vertexPairList;
339 
340  // Make sure enough hits either side
341  if (std::distance(clusHitPairVector.begin(),vertexItr) > minimumClusterSize && std::distance(vertexItr,clusHitPairVector.end()) > minimumClusterSize)
342  {
343  vertexPairList.emplace_back(Hit3DItrPair(clusHitPairVector.begin(),vertexItr));
344  vertexPairList.emplace_back(Hit3DItrPair(vertexItr,clusHitPairVector.end()));
345  }
346 
347  storeCurrentCluster = false;
348 
349  // Ok, now loop through our pairs
350  for(auto& hit3DItrPair : vertexPairList)
351  {
352  reco::ClusterParameters clusterParams;
353  reco::HitPairListPtr& hitPairListPtr = clusterParams.getHitPairListPtr();
354 
355  std::cout << indent << "+> -- building new cluster, size: " << std::distance(hit3DItrPair.first,hit3DItrPair.second) << std::endl;
356 
357  // size the container...
358  hitPairListPtr.resize(std::distance(hit3DItrPair.first,hit3DItrPair.second));
359 
360  // and copy the hits into the container
361  std::copy(hit3DItrPair.first,hit3DItrPair.second,hitPairListPtr.begin());
362 
363  // First stage of feature extraction runs here
364  m_pcaAlg.PCAAnalysis_3D(hitPairListPtr, clusterParams.getFullPCA());
365 
366  // Recover the new fullPCA
367  reco::PrincipalComponents& newFullPCA = clusterParams.getFullPCA();
368 
369  // Must have a valid pca
370  if (newFullPCA.getSvdOK())
371  {
372  std::cout << indent << "+> -- >> cluster has a valid Full PCA" << std::endl;
373 
374  // Need to check if the PCA direction has been reversed
375  Eigen::Vector3f fullPrimaryVec(fullPCA.getEigenVectors()[0][0],fullPCA.getEigenVectors()[0][1],fullPCA.getEigenVectors()[0][2]);
376  Eigen::Vector3f newPrimaryVec(newFullPCA.getEigenVectors()[0][0],newFullPCA.getEigenVectors()[0][1],newFullPCA.getEigenVectors()[0][2]);
377 
378  // If the PCA's are opposite the flip the axes
379  if (fullPrimaryVec.dot(newPrimaryVec) < 0.)
380  {
382 
383  eigenVectors.resize(3);
384 
385  for(size_t vecIdx = 0; vecIdx < 3; vecIdx++)
386  {
387  eigenVectors[vecIdx].resize(3,0.);
388 
389  eigenVectors[vecIdx][0] = -newFullPCA.getEigenVectors()[vecIdx][0];
390  eigenVectors[vecIdx][1] = -newFullPCA.getEigenVectors()[vecIdx][1];
391  eigenVectors[vecIdx][2] = -newFullPCA.getEigenVectors()[vecIdx][2];
392  }
393 
394  newFullPCA = reco::PrincipalComponents(true,
395  newFullPCA.getNumHitsUsed(),
396  newFullPCA.getEigenValues(),
397  eigenVectors,
398  newFullPCA.getAvePosition(),
399  newFullPCA.getAveHitDoca());
400  }
401 
402  // Set the skeleton PCA to make sure it has some value
403  clusterParams.getSkeletonPCA() = clusterParams.getFullPCA();
404 
405  positionItr = breakIntoTinyBits(clusterParams, positionItr, outputClusterList, level+4);
406  }
407  }
408  }
409  }
410 
411  // First question, are we done breaking?
412  if (storeCurrentCluster)
413  {
414  // I think this is where we fill out the rest of the parameters?
415  // Start by adding the 2D hits...
416  // See if we can avoid duplicates by temporarily transferring to a set
417  std::set<const reco::ClusterHit2D*> hitSet;
418 
419  // Loop through 3D hits to get a set of unique 2D hits
420  for(const auto& hit3D : clusterToBreak.getHitPairListPtr())
421  {
422  for(const auto& hit2D : hit3D->getHits())
423  {
424  if (hit2D) hitSet.insert(hit2D);
425  }
426  }
427 
428  // Now add these to the new cluster
429  for(const auto& hit2D : hitSet)
430  {
431  hit2D->setStatusBit(reco::ClusterHit2D::USED);
432  clusterToBreak.UpdateParameters(hit2D);
433  }
434 
435  std::cout << indent << "*********>>> storing new subcluster of size " << clusterToBreak.getHitPairListPtr().size() << std::endl;
436 
437  positionItr = outputClusterList.insert(positionItr,clusterToBreak);
438 
439  // The above points to the element, want the next element
440  positionItr++;
441  }
442  else if (inputPositionItr == positionItr) std::cout << indent << "***** DID NOT STORE A CLUSTER *****" << std::endl;
443 
444  return positionItr;
445 }
446 
447 void VoronoiPathFinder::buildConvexHull(reco::ClusterParameters& clusterParameters, int level) const
448 {
449  // set an indention string
450  std::string minuses(level/2, '-');
451  std::string indent(level/2, ' ');
452 
453  indent += minuses;
454 
455  // The plan is to build the enclosing 2D polygon around the points in the PCA plane of most spread for this cluster
456  // To do so we need to start by building a list of 2D projections onto the plane of most spread...
457  reco::PrincipalComponents& pca = clusterParameters.getFullPCA();
458 
459  // Recover the parameters from the Principal Components Analysis that we need to project and accumulate
460  Eigen::Vector3f pcaCenter(pca.getAvePosition()[0],pca.getAvePosition()[1],pca.getAvePosition()[2]);
461  Eigen::Vector3f planeVec0(pca.getEigenVectors()[0][0],pca.getEigenVectors()[0][1],pca.getEigenVectors()[0][2]);
462  Eigen::Vector3f planeVec1(pca.getEigenVectors()[1][0],pca.getEigenVectors()[1][1],pca.getEigenVectors()[1][2]);
463  Eigen::Vector3f pcaPlaneNrml(pca.getEigenVectors()[2][0],pca.getEigenVectors()[2][1],pca.getEigenVectors()[2][2]);
464 
465  // Let's get the rotation matrix from the standard coordinate system to the PCA system.
466  Eigen::Matrix3f rotationMatrix;
467 
468  rotationMatrix << planeVec0(0), planeVec0(1), planeVec0(2),
469  planeVec1(0), planeVec1(1), planeVec1(2),
470  pcaPlaneNrml(0), pcaPlaneNrml(1), pcaPlaneNrml(2);
471 
472  //dcel2d::PointList pointList;
473  using Point = std::tuple<float,float,const reco::ClusterHit3D*>;
474  using PointList = std::list<Point>;
475  PointList pointList;
476 
477  // Loop through hits and do projection to plane
478  for(const auto& hit3D : clusterParameters.getHitPairListPtr())
479  {
480  Eigen::Vector3f pcaToHitVec(hit3D->getPosition()[0] - pcaCenter(0),
481  hit3D->getPosition()[1] - pcaCenter(1),
482  hit3D->getPosition()[2] - pcaCenter(2));
483  Eigen::Vector3f pcaToHit = rotationMatrix * pcaToHitVec;
484 
485  pointList.emplace_back(dcel2d::Point(pcaToHit(0),pcaToHit(1),hit3D));
486  }
487 
488  // Sort the point vec by increasing x, then increase y
489  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);});
490 
491  // containers for finding the "best" hull...
492  std::vector<ConvexHull> convexHullVec;
493  std::vector<PointList> rejectedListVec;
494  bool increaseDepth(pointList.size() > 5);
495  float lastArea(std::numeric_limits<float>::max());
496 
497  while(increaseDepth)
498  {
499  // Get another convexHull container
500  convexHullVec.push_back(ConvexHull(pointList));
501  rejectedListVec.push_back(PointList());
502 
503  const ConvexHull& convexHull = convexHullVec.back();
504  PointList& rejectedList = rejectedListVec.back();
505  const PointList& convexHullPoints = convexHull.getConvexHull();
506 
507  increaseDepth = false;
508 
509  if (convexHull.getConvexHullArea() > 0.)
510  {
511  std::cout << indent << "-> built convex hull, 3D hits: " << pointList.size() << " with " << convexHullPoints.size() << " vertices" << ", area: " << convexHull.getConvexHullArea() << std::endl;
512  std::cout << indent << "-> -Points:";
513  for(const auto& point : convexHullPoints)
514  std::cout << " (" << std::get<0>(point) << "," << std::get<1>(point) << ")";
515  std::cout << std::endl;
516 
517  if (convexHullVec.size() < 2 || convexHull.getConvexHullArea() < 0.8 * lastArea)
518  {
519  for(auto& point : convexHullPoints)
520  {
521  pointList.remove(point);
522  rejectedList.emplace_back(point);
523  }
524  lastArea = convexHull.getConvexHullArea();
525 // increaseDepth = true;
526  }
527  }
528  }
529 
530  // do we have a valid convexHull?
531  while(!convexHullVec.empty() && convexHullVec.back().getConvexHullArea() < 0.5)
532  {
533  convexHullVec.pop_back();
534  rejectedListVec.pop_back();
535  }
536 
537  // If we found the convex hull then build edges around the region
538  if (!convexHullVec.empty())
539  {
540  size_t nRejectedTotal(0);
541  reco::HitPairListPtr hitPairListPtr = clusterParameters.getHitPairListPtr();
542 
543  for(const auto& rejectedList : rejectedListVec)
544  {
545  nRejectedTotal += rejectedList.size();
546 
547  for(const auto& rejectedPoint : rejectedList)
548  {
549  std::cout << indent << "-> -- Point is " << convexHullVec.back().findNearestDistance(rejectedPoint) << " from nearest edge" << std::endl;
550 
551  if (convexHullVec.back().findNearestDistance(rejectedPoint) > 0.5)
552  hitPairListPtr.remove(std::get<2>(rejectedPoint));
553  }
554  }
555 
556  std::cout << indent << "-> Removed " << nRejectedTotal << " leaving " << pointList.size() << "/" << hitPairListPtr.size() << " points" << std::endl;
557 
558  // Now add "edges" to the cluster to describe the convex hull (for the display)
559  reco::Hit3DToEdgeMap& edgeMap = clusterParameters.getHit3DToEdgeMap();
560  reco::EdgeList& bestEdgeList = clusterParameters.getBestEdgeList();
561 
562  Point lastPoint = convexHullVec.back().getConvexHull().front();
563 
564  for(auto& curPoint : convexHullVec.back().getConvexHull())
565  {
566  if (curPoint == lastPoint) continue;
567 
568  const reco::ClusterHit3D* lastPoint3D = std::get<2>(lastPoint);
569  const reco::ClusterHit3D* curPoint3D = std::get<2>(curPoint);
570 
571  float distBetweenPoints = (curPoint3D->getPosition()[0] - lastPoint3D->getPosition()[0]) * (curPoint3D->getPosition()[0] - lastPoint3D->getPosition()[0])
572  + (curPoint3D->getPosition()[1] - lastPoint3D->getPosition()[1]) * (curPoint3D->getPosition()[1] - lastPoint3D->getPosition()[1])
573  + (curPoint3D->getPosition()[2] - lastPoint3D->getPosition()[2]) * (curPoint3D->getPosition()[2] - lastPoint3D->getPosition()[2]);
574 
575  distBetweenPoints = std::sqrt(distBetweenPoints);
576 
577  reco::EdgeTuple edge(lastPoint3D,curPoint3D,distBetweenPoints);
578 
579  edgeMap[lastPoint3D].push_back(edge);
580  edgeMap[curPoint3D].push_back(edge);
581  bestEdgeList.emplace_back(edge);
582 
583  lastPoint = curPoint;
584  }
585  }
586 
587  return;
588 }
589 
590 
591 void VoronoiPathFinder::buildVoronoiDiagram(reco::ClusterParameters& clusterParameters, int level) const
592 {
593  // The plan is to build the enclosing 2D polygon around the points in the PCA plane of most spread for this cluster
594  // To do so we need to start by building a list of 2D projections onto the plane of most spread...
595  reco::PrincipalComponents& pca = clusterParameters.getFullPCA();
596 
597  // Recover the parameters from the Principal Components Analysis that we need to project and accumulate
598  Eigen::Vector3f pcaCenter(pca.getAvePosition()[0],pca.getAvePosition()[1],pca.getAvePosition()[2]);
599  Eigen::Vector3f planeVec0(pca.getEigenVectors()[0][0],pca.getEigenVectors()[0][1],pca.getEigenVectors()[0][2]);
600  Eigen::Vector3f planeVec1(pca.getEigenVectors()[1][0],pca.getEigenVectors()[1][1],pca.getEigenVectors()[1][2]);
601  Eigen::Vector3f pcaPlaneNrml(pca.getEigenVectors()[2][0],pca.getEigenVectors()[2][1],pca.getEigenVectors()[2][2]);
602 
603  // Leave the following code bits as an example of a more complicated way to do what we are trying to do here
604  // (but in case I want to use quaternions again!)
605  //
606  //Eigen::Vector3f unitVecX(1.,0.,0.);
607  //
608  //Eigen::Quaternionf rotationMatrix = Eigen::Quaternionf::FromTwoVectors(planeVec0,unitVecX);
609  //
610  //Eigen::Matrix3f Rmatrix = rotationMatrix.toRotationMatrix();
611  //Eigen::Matrix3f RInvMatrix = rotationMatrix.inverse().toRotationMatrix();
612 
613  // Let's get the rotation matrix from the standard coordinate system to the PCA system.
614  Eigen::Matrix3f rotationMatrix;
615 
616  rotationMatrix << planeVec0(0), planeVec0(1), planeVec0(2),
617  planeVec1(0), planeVec1(1), planeVec1(2),
618  pcaPlaneNrml(0), pcaPlaneNrml(1), pcaPlaneNrml(2);
619 
620  dcel2d::PointList pointList;
621 
622  // Loop through hits and do projection to plane
623  for(const auto& hit3D : clusterParameters.getHitPairListPtr())
624  {
625  Eigen::Vector3f pcaToHitVec(hit3D->getPosition()[0] - pcaCenter(0),
626  hit3D->getPosition()[1] - pcaCenter(1),
627  hit3D->getPosition()[2] - pcaCenter(2));
628  Eigen::Vector3f pcaToHit = rotationMatrix * pcaToHitVec;
629 
630  pointList.emplace_back(dcel2d::Point(pcaToHit(0),pcaToHit(1),hit3D));
631  }
632 
633  // Sort the point vec by increasing x, then increase y
634  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);});
635 
636  // Set up the voronoi diagram builder
637  voronoi2d::VoronoiDiagram voronoiDiagram(clusterParameters.getHalfEdgeList(),clusterParameters.getVertexList(),clusterParameters.getFaceList());
638 
639  // And make the diagram
640  voronoiDiagram.buildVoronoiDiagram(pointList);
641 
642  // Recover the voronoi diagram vertex list and the container to store them in
643  dcel2d::VertexList& vertexList = clusterParameters.getVertexList();
644 
645  // Now get the inverse of the rotation matrix so we can get the vertex positions,
646  // which lie in the plane of the two largest PCA axes, in the standard coordinate system
647  Eigen::Matrix3f rotationMatrixInv = rotationMatrix.inverse();
648 
649  // Translate and fill
650  for(auto& vertex : vertexList)
651  {
652  Eigen::Vector3f coords = rotationMatrixInv * vertex.getCoords();
653 
654  coords += pcaCenter;
655 
656  vertex.setCoords(coords);
657  }
658 
659  // Now do the Convex Hull
660  // Now add "edges" to the cluster to describe the convex hull (for the display)
661  reco::Hit3DToEdgeMap& edgeMap = clusterParameters.getHit3DToEdgeMap();
662  reco::EdgeList& bestEdgeList = clusterParameters.getBestEdgeList();
663 
664 // const dcel2d::PointList& edgePoints = voronoiDiagram.getConvexHull();
665 // PointList localList;
666 //
667 // for(const auto& edgePoint : edgePoints) localList.emplace_back(std::get<0>(edgePoint),std::get<1>(edgePoint),std::get<2>(edgePoint));
668 //
669 // // Sort the point vec by increasing x, then increase y
670 // 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);});
671 //
672 // // Why do I need to do this?
673 // ConvexHull convexHull(localList);
674 //
675 // Point lastPoint = convexHull.getConvexHull().front();
676  dcel2d::Point lastPoint = voronoiDiagram.getConvexHull().front();
677 
678 // std::cout << "@@@@>> Build convex hull, voronoi handed " << edgePoints.size() << " points, convexHull cut to " << convexHull.getConvexHull().size() << std::endl;
679 
680 // for(auto& curPoint : convexHull.getConvexHull())
681  for(auto& curPoint : voronoiDiagram.getConvexHull())
682  {
683  if (curPoint == lastPoint) continue;
684 
685  const reco::ClusterHit3D* lastPoint3D = std::get<2>(lastPoint);
686  const reco::ClusterHit3D* curPoint3D = std::get<2>(curPoint);
687 
688  float distBetweenPoints = (curPoint3D->getPosition()[0] - lastPoint3D->getPosition()[0]) * (curPoint3D->getPosition()[0] - lastPoint3D->getPosition()[0])
689  + (curPoint3D->getPosition()[1] - lastPoint3D->getPosition()[1]) * (curPoint3D->getPosition()[1] - lastPoint3D->getPosition()[1])
690  + (curPoint3D->getPosition()[2] - lastPoint3D->getPosition()[2]) * (curPoint3D->getPosition()[2] - lastPoint3D->getPosition()[2]);
691 
692  distBetweenPoints = std::sqrt(distBetweenPoints);
693 
694  reco::EdgeTuple edge(lastPoint3D,curPoint3D,distBetweenPoints);
695 
696  edgeMap[lastPoint3D].push_back(edge);
697  edgeMap[curPoint3D].push_back(edge);
698  bestEdgeList.emplace_back(edge);
699 
700  lastPoint = curPoint;
701  }
702 
703  std::cout << "****> vertexList containted " << vertexList.size() << " vertices" << std::endl;
704 
705  return;
706 }
707 
708 
709 float VoronoiPathFinder::closestApproach(const Eigen::Vector3f& P0,
710  const Eigen::Vector3f& u0,
711  const Eigen::Vector3f& P1,
712  const Eigen::Vector3f& u1,
713  Eigen::Vector3f& poca0,
714  Eigen::Vector3f& poca1) const
715 {
716  // Technique is to compute the arclength to each point of closest approach
717  Eigen::Vector3f w0 = P0 - P1;
718  float a(1.);
719  float b(u0.dot(u1));
720  float c(1.);
721  float d(u0.dot(w0));
722  float e(u1.dot(w0));
723  float den(a * c - b * b);
724 
725  float arcLen0 = (b * e - c * d) / den;
726  float arcLen1 = (a * e - b * d) / den;
727 
728  poca0 = P0 + arcLen0 * u0;
729  poca1 = P1 + arcLen1 * u1;
730 
731  return (poca0 - poca1).norm();
732 }
733 
735 } // 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
float getTimeToExecute() const override
If monitoring, recover the time to execute a particular function.
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.
intermediate_table::iterator iterator
Float_t den
Definition: plot.C:37
Declaration of signal hit object.
float closestApproach(const Eigen::Vector3f &, const Eigen::Vector3f &, const Eigen::Vector3f &, const Eigen::Vector3f &, Eigen::Vector3f &, Eigen::Vector3f &) const
reco::PrincipalComponents & getSkeletonPCA()
Definition: Cluster3D.h:406
reco::EdgeList & getBestEdgeList()
Definition: Cluster3D.h:409
std::pair< Point, Point > MinMaxPoints
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.
IClusterModAlg interface class definiton.
bool m_enableMonitoring
Data members to follow.
ClusterParametersList & daughterList()
Definition: Cluster3D.h:378
Implements a ConvexHull for use in clustering.
std::list< EdgeTuple > EdgeList
Definition: Cluster3D.h:324
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
void ModifyClusters(reco::ClusterParametersList &) const override
Scan an input collection of clusters and modify those according to the specific implementing algorith...
T get(std::string const &key) const
Definition: ParameterSet.h:231
std::tuple< float, float, const reco::ClusterHit3D * > Point
std::string indent(std::size_t const i)
void buildConvexHull(reco::ClusterParameters &clusterParameters, int level=0) const
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
const float * getPosition() const
Definition: Cluster3D.h:145
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.
constexpr auto const & left(const_AssnsIter< L, R, D, Dir > const &a, const_AssnsIter< L, R, D, Dir > const &b)
Definition: AssnsIter.h:104
std::unique_ptr< lar_cluster3d::IClusterAlg > m_clusterAlg
Algorithm to do 3D space point clustering.
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.
size_t m_minTinyClusterSize
Minimum size for a "tiny" cluster.
void buildVoronoiDiagram(reco::ClusterParameters &clusterParameters, int level=0) const
std::list< Point > PointList
Definition: DCEL.h:35
std::unordered_map< const reco::ClusterHit3D *, reco::EdgeList > Hit3DToEdgeMap
Definition: Cluster3D.h:326
std::pair< MinMaxPoints, MinMaxPoints > MinMaxPointPair
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
VoronoiPathFinder(const fhicl::ParameterSet &)
Constructor.
const EigenVectors & getEigenVectors() const
Definition: Cluster3D.h:227
dcel2d::VertexList & getVertexList()
Definition: Cluster3D.h:411
void configure(fhicl::ParameterSet const &pset) override
vertex reconstruction