All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
VoronoiPathFinder_tool.cc
Go to the documentation of this file.
1 /**
2  * @file VoronoiPathFinder_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 "art_root_io/TFileDirectory.h"
12 #include "cetlib/cpu_timer.h"
13 #include "fhiclcpp/ParameterSet.h"
14 #include "messagefacility/MessageLogger/MessageLogger.h"
15 
19 
20 // LArSoft includes
23 
24 // Eigen
25 #include <Eigen/Core>
26 
27 // Root histograms
28 #include "TH1F.h"
29 
30 // std includes
31 #include <string>
32 #include <iostream>
33 #include <memory>
34 
35 //------------------------------------------------------------------------------------------------------------------------------------------
36 // implementation follows
37 
38 namespace lar_cluster3d {
39 
40 class VoronoiPathFinder : virtual public IClusterModAlg
41 {
42 public:
43  /**
44  * @brief Constructor
45  *
46  * @param pset
47  */
48  explicit VoronoiPathFinder(const fhicl::ParameterSet&);
49 
50  /**
51  * @brief Destructor
52  */
54 
55  void configure(fhicl::ParameterSet const &pset) override;
56 
57  /**
58  * @brief Interface for initializing histograms if they are desired
59  * Note that the idea is to put hisgtograms in a subfolder
60  *
61  * @param TFileDirectory - the folder to store the hists in
62  */
63  void initializeHistograms(art::TFileDirectory&) override;
64 
65  /**
66  * @brief Scan an input collection of clusters and modify those according
67  * to the specific implementing algorithm
68  *
69  * @param clusterParametersList A list of cluster objects (parameters from associated hits)
70  */
71  void ModifyClusters(reco::ClusterParametersList&) const override;
72 
73  /**
74  * @brief If monitoring, recover the time to execute a particular function
75  */
76  float getTimeToExecute() const override {return fTimeToProcess;}
77 
78 private:
79 
80  /**
81  * @brief Use PCA to try to find path in cluster
82  *
83  * @param clusterParameters The given cluster parameters object to try to split
84  * @param clusterParametersList The list of clusters
85  */
86  reco::ClusterParametersList::iterator breakIntoTinyBits(reco::ClusterParameters& cluster,
88  reco::ClusterParametersList::iterator positionItr,
89  reco::ClusterParametersList& outputClusterList,
90  int level = 0) const;
91 
92  /**
93  * @brief Use PCA to try to find path in cluster
94  *
95  * @param clusterParameters The given cluster parameters object to try to split
96  * @param clusterParametersList The list of clusters
97  */
98  reco::ClusterParametersList::iterator subDivideCluster(reco::ClusterParameters& cluster,
100  reco::ClusterParametersList::iterator positionItr,
101  reco::ClusterParametersList& outputClusterList,
102  int level = 0) const;
103 
104  bool makeCandidateCluster(Eigen::Vector3f&,
106  reco::HitPairListPtr::iterator,
107  reco::HitPairListPtr::iterator,
108  int) const;
109 
110  float closestApproach(const Eigen::Vector3f&,
111  const Eigen::Vector3f&,
112  const Eigen::Vector3f&,
113  const Eigen::Vector3f&,
114  Eigen::Vector3f&,
115  Eigen::Vector3f&) const;
116 
117  using Point = std::tuple<float,float,const reco::ClusterHit3D*>;
118  using PointList = std::list<Point>;
119  using MinMaxPoints = std::pair<Point,Point>;
120  using MinMaxPointPair = std::pair<MinMaxPoints,MinMaxPoints>;
121 
122  void buildConvexHull(reco::ClusterParameters& clusterParameters, int level = 0) const;
123  void buildVoronoiDiagram(reco::ClusterParameters& clusterParameters, int level = 0) const;
124 
126  /**
127  * @brief FHICL parameters
128  */
129  bool fEnableMonitoring; ///<
130  size_t fMinTinyClusterSize; ///< Minimum size for a "tiny" cluster
131  mutable float fTimeToProcess; ///<
132 
133  /**
134  * @brief Histogram definitions
135  */
137 
144 
155 
156  /**
157  * @brief Tools
158  */
159  std::unique_ptr<lar_cluster3d::IClusterAlg> fClusterAlg; ///< Algorithm to do 3D space point clustering
160  PrincipalComponentsAlg fPCAAlg; // For running Principal Components Analysis
161 };
162 
163 VoronoiPathFinder::VoronoiPathFinder(fhicl::ParameterSet const &pset) :
164  fPCAAlg(pset.get<fhicl::ParameterSet>("PrincipalComponentsAlg"))
165 {
166  this->configure(pset);
167 }
168 
169 //------------------------------------------------------------------------------------------------------------------------------------------
170 
172 {
173 }
174 
175 //------------------------------------------------------------------------------------------------------------------------------------------
176 
177 void VoronoiPathFinder::configure(fhicl::ParameterSet const &pset)
178 {
179  fEnableMonitoring = pset.get<bool> ("EnableMonitoring", true );
180  fMinTinyClusterSize = pset.get<size_t>("MinTinyClusterSize",40);
181  fClusterAlg = art::make_tool<lar_cluster3d::IClusterAlg>(pset.get<fhicl::ParameterSet>("ClusterAlg"));
182 
183  fTimeToProcess = 0.;
184 
185  return;
186 }
187 
188 void VoronoiPathFinder::initializeHistograms(art::TFileDirectory& histDir)
189 {
190  // It is assumed that the input TFileDirectory has been set up to group histograms into a common
191  // folder at the calling routine's level. Here we create one more level of indirection to keep
192  // histograms made by this tool separate.
193  fFillHistograms = true;
194 
195  std::string dirName = "VoronoiPath";
196 
197  art::TFileDirectory dir = histDir.mkdir(dirName.c_str());
198 
199  // Divide into two sets of hists... those for the overall cluster and
200  // those for the subclusters
201  fTopNum3DHits = dir.make<TH1F>("TopNum3DHits", "Number 3D Hits", 200, 0., 200.);
202  fTopNumEdges = dir.make<TH1F>("TopNumEdges", "Number Edges", 200, 0., 200.);
203  fTopEigen21Ratio = dir.make<TH1F>("TopEigen21Rat", "Eigen 2/1 Ratio", 100, 0., 1.);
204  fTopEigen20Ratio = dir.make<TH1F>("TopEigen20Rat", "Eigen 2/0 Ratio", 100, 0., 1.);
205  fTopEigen10Ratio = dir.make<TH1F>("TopEigen10Rat", "Eigen 1/0 Ratio", 100, 0., 1.);
206  fTopPrimaryLength = dir.make<TH1F>("TopPrimaryLen", "Primary Length", 200, 0., 200.);
207 
208  fSubNum3DHits = dir.make<TH1F>("SubNum3DHits", "Number 3D Hits", 200, 0., 200.);
209  fSubNumEdges = dir.make<TH1F>("SubNumEdges", "Number Edges", 200, 0., 200.);
210  fSubEigen21Ratio = dir.make<TH1F>("SubEigen21Rat", "Eigen 2/1 Ratio", 100, 0., 1.);
211  fSubEigen20Ratio = dir.make<TH1F>("SubEigen20Rat", "Eigen 2/0 Ratio", 100, 0., 1.);
212  fSubEigen10Ratio = dir.make<TH1F>("SubEigen10Rat", "Eigen 1/0 Ratio", 100, 0., 1.);
213  fSubPrimaryLength = dir.make<TH1F>("SubPrimaryLen", "Primary Length", 200, 0., 200.);
214  fSubCosToPrevPCA = dir.make<TH1F>("SubCosToPrev", "Cos(theta)", 101, 0., 1.01);
215  fSubCosExtToPCA = dir.make<TH1F>("SubCosExtPCA", "Cos(theta)", 102, -1.01, 1.01);
216  fSubMaxDefect = dir.make<TH1F>("SubMaxDefect", "Max Defect", 100, 0., 50.);
217  fSubUsedDefect = dir.make<TH1F>("SubUsedDefect", "Used Defect", 100, 0., 50.);
218 
219  return;
220 }
221 
223 {
224  /**
225  * @brief Top level interface for algorithm to consider pairs of clusters from the input
226  * list and determine if they are consistent with each other and, therefore, should
227  * be merged. This is done by looking at the PCA for each cluster and looking at the
228  * projection of the primary axis along the vector connecting their centers.
229  */
230 
231  // Initial clustering is done, now trim the list and get output parameters
232  cet::cpu_timer theClockBuildClusters;
233 
234  // Start clocks if requested
235  if (fEnableMonitoring) theClockBuildClusters.start();
236 
237  int countClusters(0);
238 
239  // This is the loop over candidate 3D clusters
240  // Note that it might be that the list of candidate clusters is modified by splitting
241  // So we use the following construct to make sure we get all of them
242  reco::ClusterParametersList::iterator clusterParametersListItr = clusterParametersList.begin();
243 
244  while(clusterParametersListItr != clusterParametersList.end())
245  {
246  // Dereference to get the cluster paramters
247  reco::ClusterParameters& clusterParameters = *clusterParametersListItr;
248 
249  std::cout << "**> Looking at Cluster " << countClusters++ << ", # hits: " << clusterParameters.getHitPairListPtr().size() << std::endl;
250 
251  // It turns out that computing the convex hull surrounding the points in the 2D projection onto the
252  // plane of largest spread in the PCA is a good way to break up the cluster... and we do it here since
253  // we (currently) want this to be part of the standard output
254  buildVoronoiDiagram(clusterParameters);
255 
256  // Make sure our cluster has enough hits...
257  if (clusterParameters.getHitPairListPtr().size() > fMinTinyClusterSize)
258  {
259  // Get an interim cluster list
260  reco::ClusterParametersList reclusteredParameters;
261 
262  // Call the main workhorse algorithm for building the local version of candidate 3D clusters
263  //******** Remind me why we need to call this at this point when the same hits will be used? ********
264  //fClusterAlg->Cluster3DHits(clusterParameters.getHitPairListPtr(), reclusteredParameters);
265  reclusteredParameters.push_back(clusterParameters);
266 
267  std::cout << ">>>>>>>>>>> Reclustered to " << reclusteredParameters.size() << " Clusters <<<<<<<<<<<<<<<" << std::endl;
268 
269  // Only process non-empty results
270  if (!reclusteredParameters.empty())
271  {
272  // Loop over the reclustered set
273  for (auto& cluster : reclusteredParameters)
274  {
275  std::cout << "****> Calling breakIntoTinyBits with " << cluster.getHitPairListPtr().size() << " hits" << std::endl;
276 
277  // It turns out that computing the convex hull surrounding the points in the 2D projection onto the
278  // plane of largest spread in the PCA is a good way to break up the cluster... and we do it here since
279  // we (currently) want this to be part of the standard output
281 
282  // Break our cluster into smaller elements...
283  subDivideCluster(cluster, cluster.getFullPCA(), cluster.daughterList().end(), cluster.daughterList(), 4);
284 
285  std::cout << "****> Broke Cluster with " << cluster.getHitPairListPtr().size() << " into " << cluster.daughterList().size() << " sub clusters";
286  for(auto& clus : cluster.daughterList()) std::cout << ", " << clus.getHitPairListPtr().size();
287  std::cout << std::endl;
288 
289  // Add the daughters to the cluster
290  clusterParameters.daughterList().insert(clusterParameters.daughterList().end(),cluster);
291 
292  // If filling histograms we do the main cluster here
293  if (fFillHistograms)
294  {
295  reco::PrincipalComponents& fullPCA = cluster.getFullPCA();
296  std::vector<double> eigenValVec = {3. * std::sqrt(fullPCA.getEigenValues()[0]),
297  3. * std::sqrt(fullPCA.getEigenValues()[1]),
298  3. * std::sqrt(fullPCA.getEigenValues()[2])};
299  double eigen2To1Ratio = eigenValVec[0] / eigenValVec[1];
300  double eigen1To0Ratio = eigenValVec[1] / eigenValVec[2];
301  double eigen2To0Ratio = eigenValVec[0] / eigenValVec[2];
302  int num3DHits = cluster.getHitPairListPtr().size();
303  int numEdges = cluster.getBestEdgeList().size();
304 
305  fTopNum3DHits->Fill(std::min(num3DHits,199), 1.);
306  fTopNumEdges->Fill(std::min(numEdges,199), 1.);
307  fTopEigen21Ratio->Fill(eigen2To1Ratio, 1.);
308  fTopEigen20Ratio->Fill(eigen2To0Ratio, 1.);
309  fTopEigen10Ratio->Fill(eigen1To0Ratio, 1.);
310  fTopPrimaryLength->Fill(std::min(eigenValVec[0],199.), 1.);
311  }
312  }
313  }
314  }
315 
316  // Go to next cluster parameters object
317  clusterParametersListItr++;
318  }
319 
320  if (fEnableMonitoring)
321  {
322  theClockBuildClusters.stop();
323 
324  fTimeToProcess = theClockBuildClusters.accumulated_real_time();
325  }
326 
327  mf::LogDebug("Cluster3D") << ">>>>> Cluster Path finding done" << std::endl;
328 
329  return;
330 }
331 
332 reco::ClusterParametersList::iterator VoronoiPathFinder::breakIntoTinyBits(reco::ClusterParameters& clusterToBreak,
333  reco::PrincipalComponents& lastPCA,
334  reco::ClusterParametersList::iterator positionItr,
335  reco::ClusterParametersList& outputClusterList,
336  int level) const
337 {
338  // This needs to be a recursive routine...
339  // Idea is to take the input cluster and order 3D hits by arclength along PCA primary axis
340  // If the cluster is above the minimum number of hits then we divide into two and call ourself
341  // with the two halves. This means we form a new cluster with hits and PCA and then call ourself
342  // If the cluster is below the minimum then we can't break any more, simply add this cluster to
343  // the new output list.
344 
345  // set an indention string
346  std::string pluses(level/2, '+');
347  std::string indent(level/2, ' ');
348 
349  indent += pluses;
350 
351  reco::ClusterParametersList::iterator inputPositionItr = positionItr;
352 
353  // To make best use of this we'll also want the PCA for this cluster... so...
354  // Recover the prime ingredients
355  reco::PrincipalComponents& fullPCA = clusterToBreak.getFullPCA();
356  std::vector<double> eigenValVec = {3. * std::sqrt(fullPCA.getEigenValues()[0]),
357  3. * std::sqrt(fullPCA.getEigenValues()[1]),
358  3. * std::sqrt(fullPCA.getEigenValues()[2])};
359  Eigen::Vector3f fullPrimaryVec(fullPCA.getEigenVectors().row(2));
360  Eigen::Vector3f lastPrimaryVec(lastPCA.getEigenVectors().row(2));
361 
362  double cosNewToLast = std::abs(fullPrimaryVec.dot(lastPrimaryVec));
363  double eigen2To1Ratio = eigenValVec[2] / eigenValVec[1];
364  double eigen1To0Ratio = eigenValVec[1] / eigenValVec[2];
365  double eigen2To0Ratio = eigenValVec[0] / eigenValVec[2];
366  double eigen2And1Ave = 0.5 * (eigenValVec[1] + eigenValVec[0]);
367  double eigenAveTo0Ratio = eigen2And1Ave / eigenValVec[2];
368 
369  bool storeCurrentCluster(true);
370  int minimumClusterSize(fMinTinyClusterSize);
371 
372  std::cout << indent << ">>> breakIntoTinyBits with " << clusterToBreak.getHitPairListPtr().size() << " input hits, " << clusterToBreak.getBestEdgeList().size() << " edges, rat21: " << eigen2To1Ratio << ", rat20: " << eigen2To0Ratio << ", rat10: " << eigen1To0Ratio << ", ave0: " << eigenAveTo0Ratio << std::endl;
373  std::cout << indent << " --> eigen 0/1/2: " << eigenValVec[0] << "/" << eigenValVec[1] << "/" << eigenValVec[2] << ", cos: " << cosNewToLast << std::endl;
374 
375  // Create a rough cut intended to tell us when we've reached the land of diminishing returns
376  if (clusterToBreak.getBestEdgeList().size() > 5 && cosNewToLast > 0.25 && eigen2To1Ratio < 0.9 && eigen2To0Ratio > 0.001 &&
377  clusterToBreak.getHitPairListPtr().size() > size_t(2 * minimumClusterSize))
378  {
379  // We want to find the convex hull vertices that lie furthest from the line to/from the extreme points
380  // To find these we:
381  // 1) recover the extreme points
382  // 2) form the vector between them
383  // 3) loop through the vertices and keep track of distance to this vector
384  // 4) Sort the resulting list by furthest points and select the one we want
385  reco::ProjectedPointList::const_iterator extremePointListItr = clusterToBreak.getConvexHull().getConvexHullExtremePoints().begin();
386 
387  const reco::ClusterHit3D* firstEdgeHit = std::get<2>(*extremePointListItr++);
388  const reco::ClusterHit3D* secondEdgeHit = std::get<2>(*extremePointListItr);
389  Eigen::Vector3f edgeVec(secondEdgeHit->getPosition()[0] - firstEdgeHit->getPosition()[0],
390  secondEdgeHit->getPosition()[1] - firstEdgeHit->getPosition()[1],
391  secondEdgeHit->getPosition()[2] - firstEdgeHit->getPosition()[2]);
392  double edgeLen = edgeVec.norm();
393 
394  // normalize it
395  edgeVec.normalize();
396 
397  // Recover the list of 3D hits associated to this cluster
398  reco::HitPairListPtr& clusHitPairVector = clusterToBreak.getHitPairListPtr();
399 
400  // Calculate the doca to the PCA primary axis for each 3D hit
401  // Importantly, this also gives us the arclength along that axis to the hit
402  fPCAAlg.PCAAnalysis_calc3DDocas(clusHitPairVector, fullPCA);
403 
404  // Sort the hits along the PCA
405  clusHitPairVector.sort([](const auto& left, const auto& right){return left->getArclenToPoca() < right->getArclenToPoca();});
406 
407  // Set up container to keep track of edges
408  using DistEdgeTuple = std::tuple<float, const reco::EdgeTuple*>;
409  using DistEdgeTupleVec = std::vector<DistEdgeTuple>;
410 
411  DistEdgeTupleVec distEdgeTupleVec;
412 
413  // Now loop through all the edges and search for the furthers point
414  for(const auto& edge : clusterToBreak.getBestEdgeList())
415  {
416  const reco::ClusterHit3D* nextEdgeHit = std::get<0>(edge); // recover the first point
417 
418  // Create vector to this point from the longest edge
419  Eigen::Vector3f hitToEdgeVec(nextEdgeHit->getPosition()[0] - firstEdgeHit->getPosition()[0],
420  nextEdgeHit->getPosition()[1] - firstEdgeHit->getPosition()[1],
421  nextEdgeHit->getPosition()[2] - firstEdgeHit->getPosition()[2]);
422 
423  // Get projection
424  float hitProjection = hitToEdgeVec.dot(edgeVec);
425 
426  // Require that the point is really "opposite" the longest edge
427  if (hitProjection > 0. && hitProjection < edgeLen)
428  {
429  Eigen::Vector3f distToHitVec = hitToEdgeVec - hitProjection * edgeVec;
430  float distToHit = distToHitVec.norm();
431 
432  distEdgeTupleVec.emplace_back(distToHit,&edge);
433  }
434  }
435 
436  std::sort(distEdgeTupleVec.begin(),distEdgeTupleVec.end(),[](const auto& left,const auto& right){return std::get<0>(left) > std::get<0>(right);});
437 
438  for(const auto& distEdgeTuple : distEdgeTupleVec)
439  {
440  const reco::EdgeTuple& edgeTuple = *std::get<1>(distEdgeTuple);
441  const reco::ClusterHit3D* edgeHit = std::get<0>(edgeTuple);
442 
443  // Now find the hit identified above as furthest away
444  reco::HitPairListPtr::iterator vertexItr = std::find(clusHitPairVector.begin(),clusHitPairVector.end(),edgeHit);
445 
446  // Make sure enough hits either side, otherwise we just keep the current cluster
447  if (vertexItr == clusHitPairVector.end() || std::distance(clusHitPairVector.begin(),vertexItr) < minimumClusterSize || std::distance(vertexItr,clusHitPairVector.end()) < minimumClusterSize) continue;
448 
449  // Now we create a list of pairs of iterators to the start and end of each subcluster
450  using Hit3DItrPair = std::pair<reco::HitPairListPtr::iterator,reco::HitPairListPtr::iterator>;
451  using VertexPairList = std::list<Hit3DItrPair>;
452 
453  VertexPairList vertexPairList;
454 
455  vertexPairList.emplace_back(Hit3DItrPair(clusHitPairVector.begin(),vertexItr));
456  vertexPairList.emplace_back(Hit3DItrPair(vertexItr,clusHitPairVector.end()));
457 
458  storeCurrentCluster = false;
459 
460  // Ok, now loop through our pairs
461  for(auto& hit3DItrPair : vertexPairList)
462  {
463  reco::ClusterParameters clusterParams;
464  reco::HitPairListPtr& hitPairListPtr = clusterParams.getHitPairListPtr();
465 
466  std::cout << indent << "+> -- building new cluster, size: " << std::distance(hit3DItrPair.first,hit3DItrPair.second) << std::endl;
467 
468  // size the container...
469  hitPairListPtr.resize(std::distance(hit3DItrPair.first,hit3DItrPair.second));
470 
471  // and copy the hits into the container
472  std::copy(hit3DItrPair.first,hit3DItrPair.second,hitPairListPtr.begin());
473 
474  // First stage of feature extraction runs here
475  fPCAAlg.PCAAnalysis_3D(hitPairListPtr, clusterParams.getFullPCA());
476 
477  // Recover the new fullPCA
478  reco::PrincipalComponents& newFullPCA = clusterParams.getFullPCA();
479 
480  // Must have a valid pca
481  if (newFullPCA.getSvdOK())
482  {
483  std::cout << indent << "+> -- >> cluster has a valid Full PCA" << std::endl;
484 
485  // If the PCA's are opposite the flip the axes
486  if (fullPrimaryVec.dot(newFullPCA.getEigenVectors().row(2)) < 0.)
487  {
488  for(size_t vecIdx = 0; vecIdx < 3; vecIdx++) newFullPCA.flipAxis(vecIdx);
489  }
490 
491  // Set the skeleton PCA to make sure it has some value
492  clusterParams.getSkeletonPCA() = clusterParams.getFullPCA();
493 
494  // Be sure to compute the oonvex hull surrounding the now broken cluster
495  buildConvexHull(clusterParams, level+2);
496 
497  positionItr = breakIntoTinyBits(clusterParams, fullPCA, positionItr, outputClusterList, level+4);
498  }
499  }
500 
501  // If successful in breaking the cluster then we are done, otherwise try to loop
502  // again taking the next furthest hit
503  break;
504  }
505  }
506 
507  // First question, are we done breaking?
508  if (storeCurrentCluster)
509  {
510  // I think this is where we fill out the rest of the parameters?
511  // Start by adding the 2D hits...
512  // See if we can avoid duplicates by temporarily transferring to a set
513  std::set<const reco::ClusterHit2D*> hitSet;
514 
515  // Loop through 3D hits to get a set of unique 2D hits
516  for(const auto& hit3D : clusterToBreak.getHitPairListPtr())
517  {
518  for(const auto& hit2D : hit3D->getHits())
519  {
520  if (hit2D) hitSet.insert(hit2D);
521  }
522  }
523 
524  // Now add these to the new cluster
525  for(const auto& hit2D : hitSet)
526  {
527  hit2D->setStatusBit(reco::ClusterHit2D::USED);
528  clusterToBreak.UpdateParameters(hit2D);
529  }
530 
531  std::cout << indent << "*********>>> storing new subcluster of size " << clusterToBreak.getHitPairListPtr().size() << std::endl;
532 
533  positionItr = outputClusterList.insert(positionItr,clusterToBreak);
534 
535  // Are we filling histograms
536  if (fFillHistograms)
537  {
538  int num3DHits = clusterToBreak.getHitPairListPtr().size();
539  int numEdges = clusterToBreak.getBestEdgeList().size();
540  Eigen::Vector3f newPrimaryVec(fullPCA.getEigenVectors().row(2));
541  Eigen::Vector3f lastPrimaryVec(lastPCA.getEigenVectors().row(2));
542  float cosToLast = newPrimaryVec.dot(lastPrimaryVec);
543 
544  fSubNum3DHits->Fill(std::min(num3DHits,199), 1.);
545  fSubNumEdges->Fill(std::min(numEdges,199), 1.);
546  fSubEigen21Ratio->Fill(eigen2To1Ratio, 1.);
547  fSubEigen20Ratio->Fill(eigen2To0Ratio, 1.);
548  fSubEigen10Ratio->Fill(eigen1To0Ratio, 1.);
549  fSubCosToPrevPCA->Fill(cosToLast, 1.);
550  fSubPrimaryLength->Fill(std::min(eigenValVec[2],199.), 1.);
551  }
552 
553  // The above points to the element, want the next element
554  positionItr++;
555  }
556  else if (inputPositionItr != positionItr)
557  {
558  std::cout << indent << "***** DID NOT STORE A CLUSTER *****" << std::endl;
559  }
560 
561  return positionItr;
562 }
563 
564 reco::ClusterParametersList::iterator VoronoiPathFinder::subDivideCluster(reco::ClusterParameters& clusterToBreak,
565  reco::PrincipalComponents& lastPCA,
566  reco::ClusterParametersList::iterator positionItr,
567  reco::ClusterParametersList& outputClusterList,
568  int level) const
569 {
570  // This is a recursive routine to divide an input cluster, according to the maximum defect point of
571  // the convex hull until we reach the point of no further improvement.
572  // The assumption is that the input cluster is fully formed already, this routine then simply
573  // divides, if successful division into two pieces it then stores the results
574 
575  // No point in doing anything if we don't have enough space points
576  if (clusterToBreak.getHitPairListPtr().size() > size_t(2 * fMinTinyClusterSize))
577  {
578  // set an indention string
579  std::string pluses(level/2, '+');
580  std::string indent(level/2, ' ');
581 
582  indent += pluses;
583 
584  // To make best use of this we'll also want the PCA for this cluster... so...
585  // Recover the prime ingredients
586  reco::PrincipalComponents& fullPCA = clusterToBreak.getFullPCA();
587  std::vector<double> eigenValVec = {3. * std::sqrt(fullPCA.getEigenValues()[0]),
588  3. * std::sqrt(fullPCA.getEigenValues()[1]),
589  3. * std::sqrt(fullPCA.getEigenValues()[2])};
590  Eigen::Vector3f fullPrimaryVec(fullPCA.getEigenVectors().row(2));
591 
592  // We want to find the convex hull vertices that lie furthest from the line to/from the extreme points
593  // To find these we:
594  // 1) recover the extreme points
595  // 2) form the vector between them
596  // 3) loop through the vertices and keep track of distance to this vector
597  // 4) Sort the resulting list by furthest points and select the one we want
598  reco::ProjectedPointList::const_iterator extremePointListItr = clusterToBreak.getConvexHull().getConvexHullExtremePoints().begin();
599 
600  const reco::ClusterHit3D* firstEdgeHit = std::get<2>(*extremePointListItr++);
601  const reco::ClusterHit3D* secondEdgeHit = std::get<2>(*extremePointListItr);
602  Eigen::Vector3f edgeVec(secondEdgeHit->getPosition()[0] - firstEdgeHit->getPosition()[0],
603  secondEdgeHit->getPosition()[1] - firstEdgeHit->getPosition()[1],
604  secondEdgeHit->getPosition()[2] - firstEdgeHit->getPosition()[2]);
605  double edgeLen = edgeVec.norm();
606 
607  // normalize it
608  edgeVec.normalize();
609 
610  // Recover the list of 3D hits associated to this cluster
611  reco::HitPairListPtr& clusHitPairVector = clusterToBreak.getHitPairListPtr();
612 
613  // Calculate the doca to the PCA primary axis for each 3D hit
614  // Importantly, this also gives us the arclength along that axis to the hit
615  fPCAAlg.PCAAnalysis_calc3DDocas(clusHitPairVector, fullPCA);
616 
617  // Sort the hits along the PCA
618  clusHitPairVector.sort([](const auto& left, const auto& right){return left->getArclenToPoca() < right->getArclenToPoca();});
619 
620  // Set up container to keep track of edges
621  using DistEdgeTuple = std::tuple<float, const reco::EdgeTuple*>;
622  using DistEdgeTupleVec = std::vector<DistEdgeTuple>;
623 
624  DistEdgeTupleVec distEdgeTupleVec;
625 
626  // Now loop through all the edges and search for the furthers point
627  for(const auto& edge : clusterToBreak.getBestEdgeList())
628  {
629  const reco::ClusterHit3D* nextEdgeHit = std::get<0>(edge); // recover the first point
630 
631  // Create vector to this point from the longest edge
632  Eigen::Vector3f hitToEdgeVec(nextEdgeHit->getPosition()[0] - firstEdgeHit->getPosition()[0],
633  nextEdgeHit->getPosition()[1] - firstEdgeHit->getPosition()[1],
634  nextEdgeHit->getPosition()[2] - firstEdgeHit->getPosition()[2]);
635 
636  // Get projection
637  float hitProjection = hitToEdgeVec.dot(edgeVec);
638 
639  // Require that the point is really "opposite" the longest edge
640  if (hitProjection > 0. && hitProjection < edgeLen)
641  {
642  Eigen::Vector3f distToHitVec = hitToEdgeVec - hitProjection * edgeVec;
643  float distToHit = distToHitVec.norm();
644 
645  distEdgeTupleVec.emplace_back(distToHit,&edge);
646  }
647  }
648 
649  std::sort(distEdgeTupleVec.begin(),distEdgeTupleVec.end(),[](const auto& left,const auto& right){return std::get<0>(left) > std::get<0>(right);});
650 
651  // Get a temporary container to hol
652  reco::ClusterParametersList tempClusterParametersList;
653  float usedDefectDist(0.);
654 
655  for(const auto& distEdgeTuple : distEdgeTupleVec)
656  {
657  const reco::EdgeTuple& edgeTuple = *std::get<1>(distEdgeTuple);
658  const reco::ClusterHit3D* edgeHit = std::get<0>(edgeTuple);
659 
660  usedDefectDist = std::get<0>(distEdgeTuple);
661 
662  // Now find the hit identified above as furthest away
663  reco::HitPairListPtr::iterator vertexItr = std::find(clusHitPairVector.begin(),clusHitPairVector.end(),edgeHit);
664 
665  // Make sure enough hits either side, otherwise we just keep the current cluster
666  if (vertexItr == clusHitPairVector.end() || std::distance(clusHitPairVector.begin(),vertexItr) < int(fMinTinyClusterSize) || std::distance(vertexItr,clusHitPairVector.end()) < int(fMinTinyClusterSize)) continue;
667 
668  tempClusterParametersList.push_back(reco::ClusterParameters());
669 
670  reco::ClusterParameters& clusterParams1 = tempClusterParametersList.back();
671 
672  if (makeCandidateCluster(fullPrimaryVec, clusterParams1, clusHitPairVector.begin(), vertexItr, level))
673  {
674  tempClusterParametersList.push_back(reco::ClusterParameters());
675 
676  reco::ClusterParameters& clusterParams2 = tempClusterParametersList.back();
677 
678  if (makeCandidateCluster(fullPrimaryVec, clusterParams2, vertexItr, clusHitPairVector.end(), level)) break;
679  }
680 
681  // If here then we could not make two valid clusters and so we try again
682  tempClusterParametersList.clear();
683  }
684 
685  // Fallback in the event of still large clusters but not defect points
686  if (tempClusterParametersList.empty())
687  {
688  std::cout << indent << "===> no cluster cands, edgeLen: " << edgeLen << ", # hits: " << clusHitPairVector.size() << ", max defect: " << std::get<0>(distEdgeTupleVec.front()) << std::endl;
689 
690  usedDefectDist = 0.;
691 
692  if (edgeLen > 20.)
693  {
694  reco::HitPairListPtr::iterator vertexItr = clusHitPairVector.begin();
695 
696  std::advance(vertexItr, clusHitPairVector.size()/2);
697 
698  tempClusterParametersList.push_back(reco::ClusterParameters());
699 
700  reco::ClusterParameters& clusterParams1 = tempClusterParametersList.back();
701 
702  if (makeCandidateCluster(fullPrimaryVec, clusterParams1, clusHitPairVector.begin(), vertexItr, level))
703  {
704  tempClusterParametersList.push_back(reco::ClusterParameters());
705 
706  reco::ClusterParameters& clusterParams2 = tempClusterParametersList.back();
707 
708  if (!makeCandidateCluster(fullPrimaryVec, clusterParams2, vertexItr, clusHitPairVector.end(), level))
709  tempClusterParametersList.clear();
710  }
711  }
712  }
713 
714  // Can only end with no candidate clusters or two so don't
715  for(auto& clusterParams : tempClusterParametersList)
716  {
717  size_t curOutputClusterListSize = outputClusterList.size();
718 
719  positionItr = subDivideCluster(clusterParams, fullPCA, positionItr, outputClusterList, level+4);
720 
721  std::cout << indent << "Output cluster list prev: " << curOutputClusterListSize << ", now: " << outputClusterList.size() << std::endl;
722 
723  // If the cluster we sent in was successfully broken then the position iterator will be shifted
724  // This means we don't want to restore the current cluster here
725  if (curOutputClusterListSize < outputClusterList.size()) continue;
726 
727  // I think this is where we fill out the rest of the parameters?
728  // Start by adding the 2D hits...
729  // See if we can avoid duplicates by temporarily transferring to a set
730  std::set<const reco::ClusterHit2D*> hitSet;
731 
732  // Loop through 3D hits to get a set of unique 2D hits
733  for(const auto& hit3D : clusterParams.getHitPairListPtr())
734  {
735  for(const auto& hit2D : hit3D->getHits())
736  {
737  if (hit2D) hitSet.insert(hit2D);
738  }
739  }
740 
741  // Now add these to the new cluster
742  for(const auto& hit2D : hitSet)
743  {
744  hit2D->setStatusBit(reco::ClusterHit2D::USED);
745  clusterParams.UpdateParameters(hit2D);
746  }
747 
748  std::cout << indent << "*********>>> storing new subcluster of size " << clusterParams.getHitPairListPtr().size() << std::endl;
749 
750  positionItr = outputClusterList.insert(positionItr,clusterParams);
751 
752  // Are we filling histograms
753  if (fFillHistograms)
754  {
755  // Recover the new fullPCA
756  reco::PrincipalComponents& newFullPCA = clusterParams.getFullPCA();
757 
758  Eigen::Vector3f newPrimaryVec(fullPCA.getEigenVectors().row(2));
759  Eigen::Vector3f lastPrimaryVec(newFullPCA.getEigenVectors().row(2));
760 
761  int num3DHits = clusterParams.getHitPairListPtr().size();
762  int numEdges = clusterParams.getBestEdgeList().size();
763  float cosToLast = newPrimaryVec.dot(lastPrimaryVec);
764  double eigen2To1Ratio = eigenValVec[0] / eigenValVec[1];
765  double eigen1To0Ratio = eigenValVec[1] / eigenValVec[2];
766  double eigen2To0Ratio = eigenValVec[2] / eigenValVec[2];
767 
768  fSubNum3DHits->Fill(std::min(num3DHits,199), 1.);
769  fSubNumEdges->Fill(std::min(numEdges,199), 1.);
770  fSubEigen21Ratio->Fill(eigen2To1Ratio, 1.);
771  fSubEigen20Ratio->Fill(eigen2To0Ratio, 1.);
772  fSubEigen10Ratio->Fill(eigen1To0Ratio, 1.);
773  fSubCosToPrevPCA->Fill(cosToLast, 1.);
774  fSubPrimaryLength->Fill(std::min(eigenValVec[2],199.), 1.);
775  fSubCosExtToPCA->Fill(fullPrimaryVec.dot(edgeVec), 1.);
776  fSubMaxDefect->Fill(std::get<0>(distEdgeTupleVec.front()), 1.);
777  fSubUsedDefect->Fill(usedDefectDist, 1.);
778  }
779 
780  // The above points to the element, want the next element
781  positionItr++;
782  }
783  }
784 
785  return positionItr;
786 }
787 
788 bool VoronoiPathFinder::makeCandidateCluster(Eigen::Vector3f& primaryPCA,
789  reco::ClusterParameters& candCluster,
790  reco::HitPairListPtr::iterator firstHitItr,
791  reco::HitPairListPtr::iterator lastHitItr,
792  int level) const
793 {
794  std::string indent(level/2, ' ');
795 
796  reco::HitPairListPtr& hitPairListPtr = candCluster.getHitPairListPtr();
797 
798  std::cout << indent << "+> -- building new cluster, size: " << std::distance(firstHitItr,lastHitItr) << std::endl;
799 
800  // size the container...
801  hitPairListPtr.resize(std::distance(firstHitItr,lastHitItr));
802 
803  // and copy the hits into the container
804  std::copy(firstHitItr,lastHitItr,hitPairListPtr.begin());
805 
806  // First stage of feature extraction runs here
807  fPCAAlg.PCAAnalysis_3D(hitPairListPtr, candCluster.getFullPCA());
808 
809  // Recover the new fullPCA
810  reco::PrincipalComponents& newFullPCA = candCluster.getFullPCA();
811 
812  // Will we want to store this cluster?
813  bool keepThisCluster(false);
814 
815  // Must have a valid pca
816  if (newFullPCA.getSvdOK())
817  {
818  std::cout << indent << "+> -- >> cluster has a valid Full PCA" << std::endl;
819 
820  // Need to check if the PCA direction has been reversed
821  Eigen::Vector3f newPrimaryVec(newFullPCA.getEigenVectors().row(2));
822 
823  // If the PCA's are opposite the flip the axes
824  if (primaryPCA.dot(newPrimaryVec) < 0.)
825  {
826  for(size_t vecIdx = 0; vecIdx < 3; vecIdx++) newFullPCA.flipAxis(vecIdx);
827 
828  newPrimaryVec = Eigen::Vector3f(newFullPCA.getEigenVectors().row(2));
829  }
830 
831  // Set the skeleton PCA to make sure it has some value
832  candCluster.getSkeletonPCA() = candCluster.getFullPCA();
833 
834  // Be sure to compute the oonvex hull surrounding the now broken cluster
835  buildConvexHull(candCluster, level+2);
836 
837  std::vector<double> eigenValVec = {3. * std::sqrt(newFullPCA.getEigenValues()[0]),
838  3. * std::sqrt(newFullPCA.getEigenValues()[1]),
839  3. * std::sqrt(newFullPCA.getEigenValues()[2])};
840  double cosNewToLast = std::abs(primaryPCA.dot(newPrimaryVec));
841  double eigen2To1Ratio = eigenValVec[2] / eigenValVec[1];
842  double eigen1To0Ratio = eigenValVec[1] / eigenValVec[2];
843  double eigen2To0Ratio = eigenValVec[0] / eigenValVec[2];
844  double eigen2And1Ave = 0.5 * (eigenValVec[1] + eigenValVec[0]);
845  double eigenAveTo0Ratio = eigen2And1Ave / eigenValVec[2];
846 
847  std::cout << indent << ">>> subDivideClusters with " << candCluster.getHitPairListPtr().size() << " input hits, " << candCluster.getBestEdgeList().size() << " edges, rat21: " << eigen2To1Ratio << ", rat20: " << eigen2To0Ratio << ", rat10: " << eigen1To0Ratio << ", ave0: " << eigenAveTo0Ratio << std::endl;
848  std::cout << indent << " --> eigen 0/1/2: " << eigenValVec[0] << "/" << eigenValVec[1] << "/" << eigenValVec[2] << ", cos: " << cosNewToLast << std::endl;
849 
850  // Create a rough cut intended to tell us when we've reached the land of diminishing returns
851 // if (candCluster.getBestEdgeList().size() > 4 && cosNewToLast > 0.25 && eigen2To1Ratio < 0.9 && eigen2To0Ratio > 0.001)
852  if (candCluster.getBestEdgeList().size() > 4 && cosNewToLast > 0.25 && eigen2To1Ratio > 0.01 && eigen2To1Ratio < 0.99 && eigen1To0Ratio < 0.5)
853  {
854  keepThisCluster = true;
855  }
856  }
857 
858  return keepThisCluster;
859 }
860 
861 
862 void VoronoiPathFinder::buildConvexHull(reco::ClusterParameters& clusterParameters, int level) const
863 {
864  // set an indention string
865  std::string minuses(level/2, '-');
866  std::string indent(level/2, ' ');
867 
868  indent += minuses;
869 
870  // The plan is to build the enclosing 2D polygon around the points in the PCA plane of most spread for this cluster
871  // To do so we need to start by building a list of 2D projections onto the plane of most spread...
872  reco::PrincipalComponents& pca = clusterParameters.getFullPCA();
873 
874  // Recover the parameters from the Principal Components Analysis that we need to project and accumulate
875  Eigen::Vector3f pcaCenter(pca.getAvePosition()[0],pca.getAvePosition()[1],pca.getAvePosition()[2]);
876 
877  //dcel2d::PointList pointList;
878  using Point = std::tuple<float,float,const reco::ClusterHit3D*>;
879  using PointList = std::list<Point>;
880 
881  reco::ConvexHull& convexHull = clusterParameters.getConvexHull();
882  PointList pointList = convexHull.getProjectedPointList();
883 
884  // Loop through hits and do projection to plane
885  for(const auto& hit3D : clusterParameters.getHitPairListPtr())
886  {
887  Eigen::Vector3f pcaToHitVec(hit3D->getPosition()[0] - pcaCenter(0),
888  hit3D->getPosition()[1] - pcaCenter(1),
889  hit3D->getPosition()[2] - pcaCenter(2));
890  Eigen::Vector3f pcaToHit = pca.getEigenVectors() * pcaToHitVec;
891 
892  pointList.emplace_back(dcel2d::Point(pcaToHit(0),pcaToHit(1),hit3D));
893  }
894 
895  // Sort the point vec by increasing x, then increase y
896  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);});
897 
898  // containers for finding the "best" hull...
899  std::vector<ConvexHull> convexHullVec;
900  std::vector<PointList> rejectedListVec;
901  bool increaseDepth(pointList.size() > 5);
902  float lastArea(std::numeric_limits<float>::max());
903 
904  while(increaseDepth)
905  {
906  // Get another convexHull container
907  convexHullVec.push_back(ConvexHull(pointList));
908  rejectedListVec.push_back(PointList());
909 
910  const ConvexHull& convexHull = convexHullVec.back();
911  PointList& rejectedList = rejectedListVec.back();
912  const PointList& convexHullPoints = convexHull.getConvexHull();
913 
914  increaseDepth = false;
915 
916  if (convexHull.getConvexHullArea() > 0.)
917  {
918  std::cout << indent << "-> built convex hull, 3D hits: " << pointList.size() << " with " << convexHullPoints.size() << " vertices" << ", area: " << convexHull.getConvexHullArea() << std::endl;
919  std::cout << indent << "-> -Points:";
920  for(const auto& point : convexHullPoints)
921  std::cout << " (" << std::get<0>(point) << "," << std::get<1>(point) << ")";
922  std::cout << std::endl;
923 
924  if (convexHullVec.size() < 2 || convexHull.getConvexHullArea() < 0.8 * lastArea)
925  {
926  for(auto& point : convexHullPoints)
927  {
928  pointList.remove(point);
929  rejectedList.emplace_back(point);
930  }
931  lastArea = convexHull.getConvexHullArea();
932 // increaseDepth = true;
933  }
934  }
935  }
936 
937  // do we have a valid convexHull?
938  while(!convexHullVec.empty() && convexHullVec.back().getConvexHullArea() < 0.5)
939  {
940  convexHullVec.pop_back();
941  rejectedListVec.pop_back();
942  }
943 
944  // If we found the convex hull then build edges around the region
945  if (!convexHullVec.empty())
946  {
947  size_t nRejectedTotal(0);
948  reco::HitPairListPtr hitPairListPtr = clusterParameters.getHitPairListPtr();
949 
950  for(const auto& rejectedList : rejectedListVec)
951  {
952  nRejectedTotal += rejectedList.size();
953 
954  for(const auto& rejectedPoint : rejectedList)
955  {
956  std::cout << indent << "-> -- Point is " << convexHullVec.back().findNearestDistance(rejectedPoint) << " from nearest edge" << std::endl;
957 
958  if (convexHullVec.back().findNearestDistance(rejectedPoint) > 0.5)
959  hitPairListPtr.remove(std::get<2>(rejectedPoint));
960  }
961  }
962 
963  std::cout << indent << "-> Removed " << nRejectedTotal << " leaving " << pointList.size() << "/" << hitPairListPtr.size() << " points" << std::endl;
964 
965  // Now add "edges" to the cluster to describe the convex hull (for the display)
966  reco::Hit3DToEdgeMap& edgeMap = convexHull.getConvexHullEdgeMap();
967  reco::EdgeList& bestEdgeList = convexHull.getConvexHullEdgeList();
968 
969  Point lastPoint = convexHullVec.back().getConvexHull().front();
970 
971  for(auto& curPoint : convexHullVec.back().getConvexHull())
972  {
973  if (curPoint == lastPoint) continue;
974 
975  const reco::ClusterHit3D* lastPoint3D = std::get<2>(lastPoint);
976  const reco::ClusterHit3D* curPoint3D = std::get<2>(curPoint);
977 
978  float distBetweenPoints = (curPoint3D->getPosition()[0] - lastPoint3D->getPosition()[0]) * (curPoint3D->getPosition()[0] - lastPoint3D->getPosition()[0])
979  + (curPoint3D->getPosition()[1] - lastPoint3D->getPosition()[1]) * (curPoint3D->getPosition()[1] - lastPoint3D->getPosition()[1])
980  + (curPoint3D->getPosition()[2] - lastPoint3D->getPosition()[2]) * (curPoint3D->getPosition()[2] - lastPoint3D->getPosition()[2]);
981 
982  distBetweenPoints = std::sqrt(distBetweenPoints);
983 
984  reco::EdgeTuple edge(lastPoint3D,curPoint3D,distBetweenPoints);
985 
986  edgeMap[lastPoint3D].push_back(edge);
987  edgeMap[curPoint3D].push_back(edge);
988  bestEdgeList.emplace_back(edge);
989 
990  lastPoint = curPoint;
991  }
992 
993  // Store the "extreme" points
994  const ConvexHull::PointList& extremePoints = convexHullVec.back().getExtremePoints();
995  reco::ProjectedPointList& extremePointList = convexHull.getConvexHullExtremePoints();
996 
997  for(const auto& point : extremePoints) extremePointList.push_back(point);
998  }
999 
1000  return;
1001 }
1002 
1003 
1004 void VoronoiPathFinder::buildVoronoiDiagram(reco::ClusterParameters& clusterParameters, int level) const
1005 {
1006  // The plan is to build the enclosing 2D polygon around the points in the PCA plane of most spread for this cluster
1007  // To do so we need to start by building a list of 2D projections onto the plane of most spread...
1008  reco::PrincipalComponents& pca = clusterParameters.getFullPCA();
1009 
1010  // Recover the parameters from the Principal Components Analysis that we need to project and accumulate
1011  Eigen::Vector3f pcaCenter(pca.getAvePosition()[0],pca.getAvePosition()[1],pca.getAvePosition()[2]);
1012 
1013  dcel2d::PointList pointList;
1014 
1015  // Loop through hits and do projection to plane
1016  for(const auto& hit3D : clusterParameters.getHitPairListPtr())
1017  {
1018  Eigen::Vector3f pcaToHitVec(hit3D->getPosition()[0] - pcaCenter(0),
1019  hit3D->getPosition()[1] - pcaCenter(1),
1020  hit3D->getPosition()[2] - pcaCenter(2));
1021  Eigen::Vector3f pcaToHit = pca.getEigenVectors() * pcaToHitVec;
1022 
1023  pointList.emplace_back(dcel2d::Point(pcaToHit(1),pcaToHit(2),hit3D));
1024  }
1025 
1026  // Sort the point vec by increasing x, then increase y
1027  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);});
1028 
1029  std::cout << " ==> Build V diagram, sorted point list contains " << pointList.size() << " hits" << std::endl;
1030 
1031  // Set up the voronoi diagram builder
1032  voronoi2d::VoronoiDiagram voronoiDiagram(clusterParameters.getHalfEdgeList(),clusterParameters.getVertexList(),clusterParameters.getFaceList());
1033 
1034  // And make the diagram
1035  voronoiDiagram.buildVoronoiDiagram(pointList);
1036 
1037  // Recover the voronoi diagram vertex list and the container to store them in
1038  dcel2d::VertexList& vertexList = clusterParameters.getVertexList();
1039 
1040  // Now get the inverse of the rotation matrix so we can get the vertex positions,
1041  // which lie in the plane of the two largest PCA axes, in the standard coordinate system
1042  Eigen::Matrix3f rotationMatrixInv = pca.getEigenVectors().inverse();
1043 
1044  // Translate and fill
1045  for(auto& vertex : vertexList)
1046  {
1047  Eigen::Vector3f coords = rotationMatrixInv * vertex.getCoords();
1048 
1049  coords += pcaCenter;
1050 
1051  vertex.setCoords(coords);
1052  }
1053 
1054  // Now do the Convex Hull
1055  // Now add "edges" to the cluster to describe the convex hull (for the display)
1056  reco::Hit3DToEdgeMap& edgeMap = clusterParameters.getHit3DToEdgeMap();
1057  reco::EdgeList& bestEdgeList = clusterParameters.getBestEdgeList();
1058 
1059 // const dcel2d::PointList& edgePoints = voronoiDiagram.getConvexHull();
1060 // PointList localList;
1061 //
1062 // for(const auto& edgePoint : edgePoints) localList.emplace_back(std::get<0>(edgePoint),std::get<1>(edgePoint),std::get<2>(edgePoint));
1063 //
1064 // // Sort the point vec by increasing x, then increase y
1065 // 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);});
1066 //
1067 // // Why do I need to do this?
1068 // ConvexHull convexHull(localList);
1069 //
1070 // Point lastPoint = convexHull.getConvexHull().front();
1071  dcel2d::Point lastPoint = voronoiDiagram.getConvexHull().front();
1072 
1073 // std::cout << "@@@@>> Build convex hull, voronoi handed " << edgePoints.size() << " points, convexHull cut to " << convexHull.getConvexHull().size() << std::endl;
1074 
1075 // for(auto& curPoint : convexHull.getConvexHull())
1076  for(auto& curPoint : voronoiDiagram.getConvexHull())
1077  {
1078  if (curPoint == lastPoint) continue;
1079 
1080  const reco::ClusterHit3D* lastPoint3D = std::get<2>(lastPoint);
1081  const reco::ClusterHit3D* curPoint3D = std::get<2>(curPoint);
1082 
1083  float distBetweenPoints = (curPoint3D->getPosition()[0] - lastPoint3D->getPosition()[0]) * (curPoint3D->getPosition()[0] - lastPoint3D->getPosition()[0])
1084  + (curPoint3D->getPosition()[1] - lastPoint3D->getPosition()[1]) * (curPoint3D->getPosition()[1] - lastPoint3D->getPosition()[1])
1085  + (curPoint3D->getPosition()[2] - lastPoint3D->getPosition()[2]) * (curPoint3D->getPosition()[2] - lastPoint3D->getPosition()[2]);
1086 
1087  distBetweenPoints = std::sqrt(distBetweenPoints);
1088 
1089  reco::EdgeTuple edge(lastPoint3D,curPoint3D,distBetweenPoints);
1090 
1091  edgeMap[lastPoint3D].push_back(edge);
1092  edgeMap[curPoint3D].push_back(edge);
1093  bestEdgeList.emplace_back(edge);
1094 
1095  lastPoint = curPoint;
1096  }
1097 
1098  std::cout << "****> vertexList containted " << vertexList.size() << " vertices for " << clusterParameters.getHitPairListPtr().size() << " hits" << std::endl;
1099 
1100  return;
1101 }
1102 
1103 
1104 float VoronoiPathFinder::closestApproach(const Eigen::Vector3f& P0,
1105  const Eigen::Vector3f& u0,
1106  const Eigen::Vector3f& P1,
1107  const Eigen::Vector3f& u1,
1108  Eigen::Vector3f& poca0,
1109  Eigen::Vector3f& poca1) const
1110 {
1111  // Technique is to compute the arclength to each point of closest approach
1112  Eigen::Vector3f w0 = P0 - P1;
1113  float a(1.);
1114  float b(u0.dot(u1));
1115  float c(1.);
1116  float d(u0.dot(w0));
1117  float e(u1.dot(w0));
1118  float den(a * c - b * b);
1119 
1120  float arcLen0 = (b * e - c * d) / den;
1121  float arcLen1 = (a * e - b * d) / den;
1122 
1123  poca0 = P0 + arcLen0 * u0;
1124  poca1 = P1 + arcLen1 * u1;
1125 
1126  return (poca0 - poca1).norm();
1127 }
1128 
1129 float VoronoiPathFinder::findConvexHullEndPoints(const reco::EdgeList& convexHull, const reco::ClusterHit3D* first3D, const reco::ClusterHit3D* last3D) const
1130 {
1131  float largestDistance(0.);
1132 
1133  // Note that edges are vectors and that the convex hull edge list will be ordered
1134  // The idea is that the maximum distance from a given edge is to the edge just before the edge that "turns back" towards the current edge
1135  // meaning that the dot product of the two edges becomes negative.
1136  reco::EdgeList::const_iterator firstEdgeItr = convexHull.begin();
1137 
1138  while(firstEdgeItr != convexHull.end())
1139  {
1140  reco::EdgeList::const_iterator nextEdgeItr = firstEdgeItr;
1141 
1142 // Eigen::Vector2f firstEdgeVec(std::get<3>(*firstEdgeItr),std::get<);
1143 // Eigen::Vector2f lastPrimaryVec(lastPCA.getEigenVectors()[0][0],lastPCA.getEigenVectors()[0][1],lastPCA.getEigenVectors()[0][2]);
1144 // float cosToLast = newPrimaryVec.dot(lastPrimaryVec);
1145 
1146  while(++nextEdgeItr != convexHull.end())
1147  {
1148 
1149  }
1150  }
1151 
1152  return largestDistance;
1153 }
1154 
1155 DEFINE_ART_CLASS_TOOL(VoronoiPathFinder)
1156 } // 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
std::list< ProjectedPoint > ProjectedPointList
Definition: Cluster3D.h:352
walls no right
Definition: selectors.fcl:105
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:477
std::list< Point > PointList
The list of the projected points.
Definition: ConvexHull.h:34
reco::EdgeList & getBestEdgeList()
Definition: Cluster3D.h:480
std::pair< Point, Point > MinMaxPoints
reco::Hit3DToEdgeMap & getHit3DToEdgeMap()
Definition: Cluster3D.h:478
reco::HitPairListPtr & getHitPairListPtr()
Definition: Cluster3D.h:475
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
const EigenValues & getEigenValues() const
Definition: Cluster3D.h:245
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
reco::Hit3DToEdgeMap & getConvexHullEdgeMap()
Definition: Cluster3D.h:384
std::unique_ptr< lar_cluster3d::IClusterAlg > fClusterAlg
Tools.
std::tuple< double, double, const reco::ClusterHit3D * > Point
Definitions used by the VoronoiDiagram algorithm.
Definition: DCEL.h:44
reco::ProjectedPointList & getConvexHullExtremePoints()
Definition: Cluster3D.h:386
void buildVoronoiDiagram(const dcel2d::PointList &)
Given an input set of 2D points construct a 2D voronoi diagram.
Definition: Voronoi.cxx:133
Define a container for working with the convex hull.
Definition: Cluster3D.h:359
void ModifyClusters(reco::ClusterParametersList &) const override
Scan an input collection of clusters and modify those according to the specific implementing algorith...
void initializeHistograms(art::TFileDirectory &) override
Interface for initializing histograms if they are desired Note that the idea is to put hisgtograms in...
std::tuple< float, float, const reco::ClusterHit3D * > Point
const Eigen::Vector3f & getAvePosition() const
Definition: Cluster3D.h:247
void buildConvexHull(reco::ClusterParameters &clusterParameters, int level=0) const
reco::ClusterParametersList::iterator subDivideCluster(reco::ClusterParameters &cluster, reco::PrincipalComponents &lastPCA, reco::ClusterParametersList::iterator positionItr, reco::ClusterParametersList &outputClusterList, int level=0) const
Use PCA to try to find path in 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
tuple dir
Definition: dropbox.py:28
ConvexHull class definiton.
Definition: ConvexHull.h:27
This header file defines the interface to a principal components analysis designed to be used within ...
reco::ConvexHull & getConvexHull()
Definition: Cluster3D.h:481
reco::ClusterParametersList::iterator breakIntoTinyBits(reco::ClusterParameters &cluster, reco::PrincipalComponents &lastPCA, reco::ClusterParametersList::iterator positionItr, reco::ClusterParametersList &outputClusterList, int level=0) const
Use PCA to try to find path in cluster.
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 buildVoronoiDiagram(reco::ClusterParameters &clusterParameters, int level=0) const
bool makeCandidateCluster(Eigen::Vector3f &, reco::ClusterParameters &, reco::HitPairListPtr::iterator, reco::HitPairListPtr::iterator, int) const
size_t fMinTinyClusterSize
Minimum size for a &quot;tiny&quot; cluster.
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
std::pair< MinMaxPoints, MinMaxPoints > MinMaxPointPair
reco::ProjectedPointList & getProjectedPointList()
Definition: Cluster3D.h:382
reco::EdgeList & getConvexHullEdgeList()
Definition: Cluster3D.h:385
dcel2d::HalfEdgeList & getHalfEdgeList()
Definition: Cluster3D.h:484
std::list< Vertex > VertexList
Definition: DCEL.h:182
bool fFillHistograms
Histogram definitions.
BEGIN_PROLOG don t mess with this pandoraTrackGausCryoW true
BEGIN_PROLOG could also be cout
std::list< ClusterParameters > ClusterParametersList
Definition: Cluster3D.h:403
VoronoiPathFinder(const fhicl::ParameterSet &)
Constructor.
const EigenVectors & getEigenVectors() const
Definition: Cluster3D.h:246
dcel2d::VertexList & getVertexList()
Definition: Cluster3D.h:483
float getTimeToExecute() const override
If monitoring, recover the time to execute a particular function.
void configure(fhicl::ParameterSet const &pset) override
Interface for configuring the particular algorithm tool.
float findConvexHullEndPoints(const reco::EdgeList &, const reco::ClusterHit3D *, const reco::ClusterHit3D *) const