All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
ConvexHullPathFinder_tool.cc
Go to the documentation of this file.
1 /**
2  * @file ConvexHullPathFinder_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 
18 
19 // LArSoft includes
22 
23 // Eigen
24 #include <Eigen/Core>
25 
26 // Root histograms
27 #include "TH1F.h"
28 
29 // std includes
30 #include <string>
31 #include <iostream>
32 #include <memory>
33 
34 //------------------------------------------------------------------------------------------------------------------------------------------
35 // implementation follows
36 
37 namespace lar_cluster3d {
38 
39 class ConvexHullPathFinder : virtual public IClusterModAlg
40 {
41 public:
42  /**
43  * @brief Constructor
44  *
45  * @param pset
46  */
47  explicit ConvexHullPathFinder(const fhicl::ParameterSet&);
48 
49  /**
50  * @brief Destructor
51  */
53 
54  void configure(fhicl::ParameterSet const &pset) override;
55 
56  /**
57  * @brief Interface for initializing histograms if they are desired
58  * Note that the idea is to put hisgtograms in a subfolder
59  *
60  * @param TFileDirectory - the folder to store the hists in
61  */
62  void initializeHistograms(art::TFileDirectory&) override;
63 
64  /**
65  * @brief Scan an input collection of clusters and modify those according
66  * to the specific implementing algorithm
67  *
68  * @param clusterParametersList A list of cluster objects (parameters from associated hits)
69  */
70  void ModifyClusters(reco::ClusterParametersList&) const override;
71 
72  /**
73  * @brief If monitoring, recover the time to execute a particular function
74  */
75  float getTimeToExecute() const override {return fTimeToProcess;}
76 
77 private:
78 
79  /**
80  * @brief Use PCA to try to find path in cluster
81  *
82  * @param clusterParameters The given cluster parameters object to try to split
83  * @param clusterParametersList The list of clusters
84  */
85  reco::ClusterParametersList::iterator subDivideCluster(reco::ClusterParameters& cluster,
87  reco::ClusterParametersList::iterator positionItr,
88  reco::ClusterParametersList& outputClusterList,
89  int level = 0) const;
90 
91  using HitOrderTuple = std::tuple<float,float,reco::ProjectedPoint>;
92  using HitOrderTupleList = std::list<HitOrderTuple>;
93 
94  bool makeCandidateCluster(Eigen::Vector3f&,
96  reco::HitPairListPtr::iterator,
97  reco::HitPairListPtr::iterator,
98  int) const;
99 
100  bool makeCandidateCluster(Eigen::Vector3f&,
103  int) const;
104 
105  bool completeCandidateCluster(Eigen::Vector3f&, reco::ClusterParameters&, int) const;
106 
112 
113  float closestApproach(const Eigen::Vector3f&,
114  const Eigen::Vector3f&,
115  const Eigen::Vector3f&,
116  const Eigen::Vector3f&,
117  Eigen::Vector3f&,
118  Eigen::Vector3f&) const;
119 
120  using MinMaxPoints = std::pair<reco::ProjectedPoint,reco::ProjectedPoint>;
121  using MinMaxPointPair = std::pair<MinMaxPoints,MinMaxPoints>;
122 
123  void buildConvexHull(reco::ClusterParameters& clusterParameters, int level = 0) const;
124 
126 
127  using KinkTuple = std::tuple<int, reco::ConvexHullKinkTuple, HitOrderTupleList, HitOrderTupleList>;
128  using KinkTupleVec = std::vector<KinkTuple>;
129 
130  void orderHitsAlongEdge(const reco::ProjectedPointList&, const reco::ProjectedPoint&, const Eigen::Vector2f&, HitOrderTupleList&) const;
131 
133 
134  void fillConvexHullHists(reco::ClusterParameters&, bool) const;
135 
136  /**
137  * @brief FHICL parameters
138  */
139  bool fEnableMonitoring; ///<
140  size_t fMinTinyClusterSize; ///< Minimum size for a "tiny" cluster
141  float fMinGapSize; ///< Minimum gap size to break at gaps
142  float fMinEigen0To1Ratio; ///< Minimum ratio of eigen 0 to 1 to continue breaking
143  float fConvexHullKinkAngle; ///< Angle to declare a kink in convex hull calc
144  float fConvexHullMinSep; ///< Min hit separation to conisder in convex hull
145  mutable float fTimeToProcess; ///<
146 
147  /**
148  * @brief Histogram definitions
149  */
151 
161 
174 
175  /**
176  * @brief Tools
177  */
178  std::unique_ptr<lar_cluster3d::IClusterAlg> fClusterAlg; ///< Algorithm to do 3D space point clustering
179  PrincipalComponentsAlg fPCAAlg; // For running Principal Components Analysis
180 };
181 
182 ConvexHullPathFinder::ConvexHullPathFinder(fhicl::ParameterSet const &pset) :
183  fPCAAlg(pset.get<fhicl::ParameterSet>("PrincipalComponentsAlg"))
184 {
185  this->configure(pset);
186 }
187 
188 //------------------------------------------------------------------------------------------------------------------------------------------
189 
191 {
192 }
193 
194 //------------------------------------------------------------------------------------------------------------------------------------------
195 
196 void ConvexHullPathFinder::configure(fhicl::ParameterSet const &pset)
197 {
198  fEnableMonitoring = pset.get<bool> ("EnableMonitoring", true );
199  fMinTinyClusterSize = pset.get<size_t>("MinTinyClusterSize", 40 );
200  fMinGapSize = pset.get<float >("MinClusterGapSize", 2.0 );
201  fMinEigen0To1Ratio = pset.get<float >("MinEigen0To1Ratio", 10.0 );
202  fConvexHullKinkAngle = pset.get<float >("ConvexHullKinkAgle", 0.95);
203  fConvexHullMinSep = pset.get<float >("ConvexHullMinSep", 0.65);
204  fClusterAlg = art::make_tool<lar_cluster3d::IClusterAlg>(pset.get<fhicl::ParameterSet>("ClusterAlg"));
205 
206  fTimeToProcess = 0.;
207 
208  return;
209 }
210 
211 void ConvexHullPathFinder::initializeHistograms(art::TFileDirectory& histDir)
212 {
213  // It is assumed that the input TFileDirectory has been set up to group histograms into a common
214  // folder at the calling routine's level. Here we create one more level of indirection to keep
215  // histograms made by this tool separate.
216  fFillHistograms = true;
217 
218  std::string dirName = "ConvexHullPath";
219 
220  art::TFileDirectory dir = histDir.mkdir(dirName.c_str());
221 
222  // Divide into two sets of hists... those for the overall cluster and
223  // those for the subclusters
224  fTopNum3DHits = dir.make<TH1F>("TopNum3DHits", "Number 3D Hits", 200, 0., 200.);
225  fTopNumEdges = dir.make<TH1F>("TopNumEdges", "Number Edges", 200, 0., 200.);
226  fTopEigen21Ratio = dir.make<TH1F>("TopEigen21Rat", "Eigen 2/1 Ratio", 100, 0., 1.);
227  fTopEigen20Ratio = dir.make<TH1F>("TopEigen20Rat", "Eigen 2/0 Ratio", 100, 0., 1.);
228  fTopEigen10Ratio = dir.make<TH1F>("TopEigen10Rat", "Eigen 1/0 Ratio", 100, 0., 1.);
229  fTopPrimaryLength = dir.make<TH1F>("TopPrimaryLen", "Primary Length", 200, 0., 200.);
230  fTopExtremeSep = dir.make<TH1F>("TopExtremeSep", "Extreme Dist", 200, 0., 200.);
231  fTopConvexCosEdge = dir.make<TH1F>("TopConvexCos", "CH Edge Cos", 100, -1., 1.);
232  fTopConvexEdgeLen = dir.make<TH1F>("TopConvexEdge", "CH Edge Len", 200, 0., 50.);
233 
234  fSubNum3DHits = dir.make<TH1F>("SubNum3DHits", "Number 3D Hits", 200, 0., 200.);
235  fSubNumEdges = dir.make<TH1F>("SubNumEdges", "Number Edges", 200, 0., 200.);
236  fSubEigen21Ratio = dir.make<TH1F>("SubEigen21Rat", "Eigen 2/1 Ratio", 100, 0., 1.);
237  fSubEigen20Ratio = dir.make<TH1F>("SubEigen20Rat", "Eigen 2/0 Ratio", 100, 0., 1.);
238  fSubEigen10Ratio = dir.make<TH1F>("SubEigen10Rat", "Eigen 1/0 Ratio", 100, 0., 1.);
239  fSubPrimaryLength = dir.make<TH1F>("SubPrimaryLen", "Primary Length", 200, 0., 200.);
240  fSubCosToPrevPCA = dir.make<TH1F>("SubCosToPrev", "Cos(theta)", 101, 0., 1.01);
241  fSubCosExtToPCA = dir.make<TH1F>("SubCosExtPCA", "Cos(theta)", 102, -1.01, 1.01);
242  fSubConvexCosEdge = dir.make<TH1F>("SubConvexCos", "CH Edge Cos", 100, -1., 1.);
243  fSubConvexEdgeLen = dir.make<TH1F>("SubConvexEdge", "CH Edge Len", 200, 0., 50.);
244  fSubMaxDefect = dir.make<TH1F>("SubMaxDefect", "Max Defect", 100, 0., 50.);
245  fSubUsedDefect = dir.make<TH1F>("SubUsedDefect", "Used Defect", 100, 0., 50.);
246 
247  return;
248 }
249 
251 {
252  /**
253  * @brief Top level interface for algorithm to consider pairs of clusters from the input
254  * list and determine if they are consistent with each other and, therefore, should
255  * be merged. This is done by looking at the PCA for each cluster and looking at the
256  * projection of the primary axis along the vector connecting their centers.
257  */
258 
259  // Initial clustering is done, now trim the list and get output parameters
260  cet::cpu_timer theClockBuildClusters;
261 
262  // Start clocks if requested
263  if (fEnableMonitoring) theClockBuildClusters.start();
264 
265  // This is the loop over candidate 3D clusters
266  // Note that it might be that the list of candidate clusters is modified by splitting
267  // So we use the following construct to make sure we get all of them
268  reco::ClusterParametersList::iterator clusterParametersListItr = clusterParametersList.begin();
269 
270  while(clusterParametersListItr != clusterParametersList.end())
271  {
272  // Dereference to get the cluster paramters
273  reco::ClusterParameters& clusterParameters = *clusterParametersListItr;
274 
275  // It turns out that computing the convex hull surrounding the points in the 2D projection onto the
276  // plane of largest spread in the PCA is a good way to break up the cluster... and we do it here since
277  // we (currently) want this to be part of the standard output
278  buildConvexHull(clusterParameters);
279 
280  // Make sure our cluster has enough hits...
281  if (clusterParameters.getHitPairListPtr().size() > fMinTinyClusterSize)
282  {
283  // Get an interim cluster list
284  reco::ClusterParametersList reclusteredParameters;
285 
286  // Call the main workhorse algorithm for building the local version of candidate 3D clusters
287  //******** Remind me why we need to call this at this point when the same hits will be used? ********
288  //fClusterAlg->Cluster3DHits(clusterParameters.getHitPairListPtr(), reclusteredParameters);
289  reclusteredParameters.push_back(clusterParameters);
290 
291  // Only process non-empty results
292  if (!reclusteredParameters.empty())
293  {
294  // Loop over the reclustered set
295  for (auto& cluster : reclusteredParameters)
296  {
297  // It turns out that computing the convex hull surrounding the points in the 2D projection onto the
298  // plane of largest spread in the PCA is a good way to break up the cluster... and we do it here since
299  // we (currently) want this to be part of the standard output
301 
302  // Break our cluster into smaller elements...
303  subDivideCluster(cluster, cluster.getFullPCA(), cluster.daughterList().end(), cluster.daughterList(), 0);
304 
305  // Add the daughters to the cluster
306  clusterParameters.daughterList().insert(clusterParameters.daughterList().end(),cluster);
307 
308  // If filling histograms we do the main cluster here
309  if (fFillHistograms)
310  {
311  reco::PrincipalComponents& fullPCA = cluster.getFullPCA();
312  std::vector<double> eigenValVec = {3. * std::sqrt(fullPCA.getEigenValues()[0]),
313  3. * std::sqrt(fullPCA.getEigenValues()[1]),
314  3. * std::sqrt(fullPCA.getEigenValues()[2])};
315  double eigen2To1Ratio = eigenValVec[0] / eigenValVec[1];
316  double eigen1To0Ratio = eigenValVec[1] / eigenValVec[2];
317  double eigen2To0Ratio = eigenValVec[2] / eigenValVec[2];
318  int num3DHits = cluster.getHitPairListPtr().size();
319  int numEdges = cluster.getConvexHull().getConvexHullEdgeList().size();
320 
321  fTopNum3DHits->Fill(std::min(num3DHits,199), 1.);
322  fTopNumEdges->Fill(std::min(numEdges,199), 1.);
323  fTopEigen21Ratio->Fill(eigen2To1Ratio, 1.);
324  fTopEigen20Ratio->Fill(eigen2To0Ratio, 1.);
325  fTopEigen10Ratio->Fill(eigen1To0Ratio, 1.);
326  fTopPrimaryLength->Fill(std::min(eigenValVec[2],199.), 1.);
327 // fTopExtremeSep->Fill(std::min(edgeLen,199.), 1.);
328  fillConvexHullHists(clusterParameters, true);
329  }
330  }
331  }
332  }
333 
334  // Go to next cluster parameters object
335  clusterParametersListItr++;
336  }
337 
338  if (fEnableMonitoring)
339  {
340  theClockBuildClusters.stop();
341 
342  fTimeToProcess = theClockBuildClusters.accumulated_real_time();
343  }
344 
345  mf::LogDebug("Cluster3D") << ">>>>> Cluster Path finding done" << std::endl;
346 
347  return;
348 }
349 
350 reco::ClusterParametersList::iterator ConvexHullPathFinder::subDivideCluster(reco::ClusterParameters& clusterToBreak,
351  reco::PrincipalComponents& lastPCA,
352  reco::ClusterParametersList::iterator positionItr,
353  reco::ClusterParametersList& outputClusterList,
354  int level) const
355 {
356  // This is a recursive routine to divide an input cluster, according to the maximum defect point of
357  // the convex hull until we reach the point of no further improvement.
358  // The assumption is that the input cluster is fully formed already, this routine then simply
359  // divides, if successful division into two pieces it then stores the results
360 
361  // No point in doing anything if we don't have enough space points
362  if (clusterToBreak.getHitPairListPtr().size() > size_t(2 * fMinTinyClusterSize))
363  {
364  // set an indention string
365 // std::string pluses(level/2, '+');
366 // std::string indent(level/2, ' ');
367 //
368 // indent += pluses;
369 
370  // We want to find the convex hull vertices that lie furthest from the line to/from the extreme points
371  // To find these we:
372  // 1) recover the extreme points
373  // 2) form the vector between them
374  // 3) loop through the vertices and keep track of distance to this vector
375  // 4) Sort the resulting list by furthest points and select the one we want
376  reco::ConvexHull& convexHull = clusterToBreak.getConvexHull();
377 
378  reco::ProjectedPointList::const_iterator extremePointListItr = convexHull.getConvexHullExtremePoints().begin();
379 
380  const reco::ClusterHit3D* firstEdgeHit = std::get<2>(*extremePointListItr++);
381  const reco::ClusterHit3D* secondEdgeHit = std::get<2>(*extremePointListItr);
382  Eigen::Vector3f edgeVec(secondEdgeHit->getPosition()[0] - firstEdgeHit->getPosition()[0],
383  secondEdgeHit->getPosition()[1] - firstEdgeHit->getPosition()[1],
384  secondEdgeHit->getPosition()[2] - firstEdgeHit->getPosition()[2]);
385 // double edgeLen = edgeVec.norm();
386 
387  // normalize it
388  edgeVec.normalize();
389 
390  // Recover the PCA for the input cluster
391  reco::PrincipalComponents& fullPCA = clusterToBreak.getFullPCA();
392  Eigen::Vector3f fullPrimaryVec(fullPCA.getEigenVectors().row(2));
393 
394  // Calculate the doca to the PCA primary axis for each 3D hit
395  // Importantly, this also gives us the arclength along that axis to the hit
396 // fPCAAlg.PCAAnalysis_calc3DDocas(clusHitPairVector, fullPCA);
397 
398  // Sort the hits along the PCA
399 // clusHitPairVector.sort([](const auto& left, const auto& right){return left->getArclenToPoca() < right->getArclenToPoca();});
400 
401  // Get a temporary container to hol
402  reco::ClusterParametersList tempClusterParametersList;
403 
404  // Try breaking clusters by finding the "maximum defect" point.
405  // If this fails the fallback in the event of still large clusters is to split in half
406  // If starting with the top level cluster then we first try to break at the kinks
407  if (level == 0)
408  {
409  // if (!breakClusterByMaxDefect(clusterToBreak, clusHitPairVector, tempClusterParametersList, level))
410  if (!breakClusterByKinks(clusterToBreak, tempClusterParametersList, level))
411  {
412  // Look to see if we can break at a gap
413  if (!breakClusterAtBigGap(clusterToBreak, tempClusterParametersList, level))
414  {
415  // It might be that we have a large deviation in the convex hull...
416  if (!breakClusterByMaxDefect(clusterToBreak, tempClusterParametersList, level))
417  {
418  std::vector<double> eigenValVec = {3. * std::sqrt(fullPCA.getEigenValues()[0]),
419  3. * std::sqrt(fullPCA.getEigenValues()[1]),
420  3. * std::sqrt(fullPCA.getEigenValues()[2])};
421 
422  // Well, we don't want "flippers" so make sure the edge has some teeth to it
423  //if (edgeLen > 10.) breakClusterInHalf(clusterToBreak, clusHitPairVector, tempClusterParametersList, level);
424  if (eigenValVec[2] > fMinEigen0To1Ratio * eigenValVec[1]) breakClusterInHalf(clusterToBreak, tempClusterParametersList, level);
425  }
426  }
427  }
428 
429  }
430  // Otherwise, change the order
431  else
432  {
433  // if (!breakClusterByMaxDefect(clusterToBreak, clusHitPairVector, tempClusterParametersList, level))
434  if (!breakClusterAtBigGap(clusterToBreak, tempClusterParametersList, level))
435  {
436  // Look to see if we can break at a gap
437  if (!breakClusterByKinks(clusterToBreak, tempClusterParametersList, level))
438  {
439  // It might be that we have a large deviation in the convex hull...
440  if (!breakClusterByMaxDefect(clusterToBreak, tempClusterParametersList, level))
441  {
442  std::vector<double> eigenValVec = {3. * std::sqrt(fullPCA.getEigenValues()[0]),
443  3. * std::sqrt(fullPCA.getEigenValues()[1]),
444  3. * std::sqrt(fullPCA.getEigenValues()[2])};
445 
446  // Well, we don't want "flippers" so make sure the edge has some teeth to it
447  //if (edgeLen > 10.) breakClusterInHalf(clusterToBreak, clusHitPairVector, tempClusterParametersList, level);
448  if (eigenValVec[2] > fMinEigen0To1Ratio * eigenValVec[1]) breakClusterInHalf(clusterToBreak, tempClusterParametersList, level);
449  }
450  }
451  }
452 
453  }
454 
455  // Can only end with no candidate clusters or two so don't
456  for(auto& clusterParams : tempClusterParametersList)
457  {
458  size_t curOutputClusterListSize = outputClusterList.size();
459 
460  positionItr = subDivideCluster(clusterParams, fullPCA, positionItr, outputClusterList, level+4);
461 
462  // If the cluster we sent in was successfully broken then the position iterator will be shifted
463  // This means we don't want to restore the current cluster here
464  if (curOutputClusterListSize < outputClusterList.size()) continue;
465 
466  // The current cluster was not further subdivided so we store its info here
467  // I think this is where we fill out the rest of the parameters?
468  // Start by adding the 2D hits...
469  // See if we can avoid duplicates by temporarily transferring to a set
470  std::set<const reco::ClusterHit2D*> hitSet;
471 
472  // Loop through 3D hits to get a set of unique 2D hits
473  for(const auto& hit3D : clusterParams.getHitPairListPtr())
474  {
475  for(const auto& hit2D : hit3D->getHits())
476  {
477  if (hit2D) hitSet.insert(hit2D);
478  }
479  }
480 
481  // Now add these to the new cluster
482  for(const auto& hit2D : hitSet)
483  {
484  hit2D->setStatusBit(reco::ClusterHit2D::USED);
485  clusterParams.UpdateParameters(hit2D);
486  }
487 
488  positionItr = outputClusterList.insert(positionItr,clusterParams);
489 
490  // Are we filling histograms
491  if (fFillHistograms)
492  {
493  std::vector<double> eigenValVec = {3. * std::sqrt(fullPCA.getEigenValues()[0]),
494  3. * std::sqrt(fullPCA.getEigenValues()[1]),
495  3. * std::sqrt(fullPCA.getEigenValues()[2])};
496 
497  // Recover the new fullPCA
498  reco::PrincipalComponents& newFullPCA = clusterParams.getFullPCA();
499 
500  Eigen::Vector3f newPrimaryVec(fullPCA.getEigenVectors().row(2));
501  Eigen::Vector3f lastPrimaryVec(newFullPCA.getEigenVectors().row(2));
502 
503  int num3DHits = clusterParams.getHitPairListPtr().size();
504  int numEdges = clusterParams.getConvexHull().getConvexHullEdgeList().size();
505  float cosToLast = newPrimaryVec.dot(lastPrimaryVec);
506  double eigen2To1Ratio = eigenValVec[0] / eigenValVec[1];
507  double eigen1To0Ratio = eigenValVec[1] / eigenValVec[2];
508  double eigen2To0Ratio = eigenValVec[0] / eigenValVec[2];
509 
510  fSubNum3DHits->Fill(std::min(num3DHits,199), 1.);
511  fSubNumEdges->Fill(std::min(numEdges,199), 1.);
512  fSubEigen21Ratio->Fill(eigen2To1Ratio, 1.);
513  fSubEigen20Ratio->Fill(eigen2To0Ratio, 1.);
514  fSubEigen10Ratio->Fill(eigen1To0Ratio, 1.);
515  fSubCosToPrevPCA->Fill(cosToLast, 1.);
516  fSubPrimaryLength->Fill(std::min(eigenValVec[2],199.), 1.);
517  fSubCosExtToPCA->Fill(fullPrimaryVec.dot(edgeVec), 1.);
518  }
519 
520  // The above points to the element, want the next element
521  positionItr++;
522  }
523  }
524 
525  return positionItr;
526 }
527 
528 bool ConvexHullPathFinder::makeCandidateCluster(Eigen::Vector3f& primaryPCA,
529  reco::ClusterParameters& candCluster,
530  reco::HitPairListPtr::iterator firstHitItr,
531  reco::HitPairListPtr::iterator lastHitItr,
532  int level) const
533 {
534  std::string indent(level/2, ' ');
535 
536  reco::HitPairListPtr& hitPairListPtr = candCluster.getHitPairListPtr();
537 
538  // size the container...
539  hitPairListPtr.resize(std::distance(firstHitItr,lastHitItr));
540 
541  // and copy the hits into the container
542  std::copy(firstHitItr,lastHitItr,hitPairListPtr.begin());
543 
544  // Will we want to store this cluster?
545  bool keepThisCluster(false);
546 
547  // Must have a valid pca
548  if (completeCandidateCluster(primaryPCA, candCluster, level))
549  {
550  // Recover the new fullPCA
551  reco::PrincipalComponents& newFullPCA = candCluster.getFullPCA();
552 
553  // Need to check if the PCA direction has been reversed
554  Eigen::Vector3f newPrimaryVec(newFullPCA.getEigenVectors().row(2));
555 
556  std::vector<double> eigenValVec = {3. * std::sqrt(newFullPCA.getEigenValues()[0]),
557  3. * std::sqrt(newFullPCA.getEigenValues()[1]),
558  3. * std::sqrt(newFullPCA.getEigenValues()[2])};
559  double cosNewToLast = std::abs(primaryPCA.dot(newPrimaryVec));
560  double eigen2To1Ratio = eigenValVec[0] / eigenValVec[1];
561  double eigen1To0Ratio = eigenValVec[1] / eigenValVec[2];
562 
563  // Create a rough cut intended to tell us when we've reached the land of diminishing returns
564 // if (candCluster.getBestEdgeList().size() > 4 && cosNewToLast > 0.25 && eigen2To1Ratio < 0.9 && eigen2To0Ratio > 0.001)
565  if (candCluster.getConvexHull().getConvexHullEdgeList().size() > 4 && cosNewToLast > 0.25 && eigen2To1Ratio > 0.01 && eigen2To1Ratio < 0.99 && eigen1To0Ratio < 0.5)
566  {
567  keepThisCluster = true;
568  }
569  }
570 
571  return keepThisCluster;
572 }
573 
574 bool ConvexHullPathFinder::makeCandidateCluster(Eigen::Vector3f& primaryPCA,
575  reco::ClusterParameters& candCluster,
576  HitOrderTupleList& orderedList,
577  int level) const
578 {
579  std::string indent(level/2, ' ');
580 
581  reco::HitPairListPtr& hitPairListPtr = candCluster.getHitPairListPtr();
582 
583  // Fill the list the old fashioned way...
584  for(const auto& tupleVal : orderedList) hitPairListPtr.emplace_back(std::get<2>(std::get<2>(tupleVal)));
585 
586  // Will we want to store this cluster?
587  bool keepThisCluster(false);
588 
589  // Must have a valid pca
590  if (completeCandidateCluster(primaryPCA, candCluster, level))
591  {
592  // Recover the new fullPCA
593  reco::PrincipalComponents& newFullPCA = candCluster.getFullPCA();
594 
595  // Need to check if the PCA direction has been reversed
596  Eigen::Vector3f newPrimaryVec(newFullPCA.getEigenVectors().row(2));
597 
598  std::vector<double> eigenValVec = {3. * std::sqrt(newFullPCA.getEigenValues()[0]),
599  3. * std::sqrt(newFullPCA.getEigenValues()[1]),
600  3. * std::sqrt(newFullPCA.getEigenValues()[2])};
601  double cosNewToLast = std::abs(primaryPCA.dot(newPrimaryVec));
602  double eigen2To1Ratio = eigenValVec[0] / eigenValVec[1];
603  double eigen1To0Ratio = eigenValVec[1] / eigenValVec[2];
604 
605  // Create a rough cut intended to tell us when we've reached the land of diminishing returns
606  // if (candCluster.getBestEdgeList().size() > 4 && cosNewToLast > 0.25 && eigen2To1Ratio < 0.9 && eigen2To0Ratio > 0.001)
607  if (candCluster.getBestEdgeList().size() > 4 && cosNewToLast > 0.25 && eigen2To1Ratio > 0.01 && eigen2To1Ratio < 0.99 && eigen1To0Ratio < 0.5)
608  {
609  keepThisCluster = true;
610  }
611  }
612 
613  return keepThisCluster;
614 }
615 
616 bool ConvexHullPathFinder::completeCandidateCluster(Eigen::Vector3f& primaryPCA, reco::ClusterParameters& candCluster, int level) const
617 {
618  // First stage of feature extraction runs here
619  fPCAAlg.PCAAnalysis_3D(candCluster.getHitPairListPtr(), candCluster.getFullPCA());
620 
621  // Recover the new fullPCA
622  reco::PrincipalComponents& newFullPCA = candCluster.getFullPCA();
623 
624  // Will we want to store this cluster?
625  bool keepThisCluster(false);
626 
627  // Must have a valid pca
628  if (newFullPCA.getSvdOK())
629  {
630  // Need to check if the PCA direction has been reversed
631  Eigen::Vector3f newPrimaryVec(newFullPCA.getEigenVectors().row(2));
632 
633  // If the PCA's are opposite the flip the axes
634  if (primaryPCA.dot(newPrimaryVec) < 0.)
635  {
636  for(size_t vecIdx = 0; vecIdx < 3; vecIdx++) newFullPCA.flipAxis(vecIdx);
637  }
638 
639  // Set the skeleton PCA to make sure it has some value
640  candCluster.getSkeletonPCA() = candCluster.getFullPCA();
641 
642  // Be sure to compute the oonvex hull surrounding the now broken cluster
643  buildConvexHull(candCluster, level+2);
644 
645  keepThisCluster = true;
646  }
647 
648  return keepThisCluster;
649 }
650 
652 {
653  // Set up container to keep track of edges
654  using HitKinkTuple = std::tuple<int, reco::HitPairListPtr::iterator>;
655  using HitKinkTupleVec = std::vector<HitKinkTuple>;
656 
657  // Recover our hits
658  reco::HitPairListPtr& hitList = clusterToBreak.getHitPairListPtr();
659 
660  // Set up container to keep track of edges
661  HitKinkTupleVec kinkTupleVec;
662 
663  reco::ConvexHull& convexHull = clusterToBreak.getConvexHull();
664  reco::ConvexHullKinkTupleList& kinkPointList = convexHull.getConvexHullKinkPoints();
665 
666  for(auto& kink : kinkPointList)
667  {
668  const reco::ClusterHit3D* hit3D = std::get<2>(std::get<0>(kink));
669 
670  reco::HitPairListPtr::iterator kinkItr = std::find(hitList.begin(),hitList.end(),hit3D);
671 
672  if (kinkItr == hitList.end()) continue;
673 
674  int numStartToKink = std::distance(hitList.begin(),kinkItr);
675  int numKinkToEnd = std::distance(kinkItr, hitList.end());
676  int minNumHits = std::min(numStartToKink,numKinkToEnd);
677 
678  if (minNumHits > int(fMinTinyClusterSize)) kinkTupleVec.emplace_back(minNumHits,kinkItr);
679  }
680 
681  // No work if the list is empty
682  if (!kinkTupleVec.empty())
683  {
684  std::sort(kinkTupleVec.begin(),kinkTupleVec.end(),[](const auto& left,const auto& right){return std::get<0>(left) > std::get<0>(right);});
685 
686  // Recover the kink point
687  reco::HitPairListPtr::iterator kinkItr = std::get<1>(kinkTupleVec.front());
688 
689  // Set up to split the input cluster
690  outputClusterList.push_back(reco::ClusterParameters());
691 
692  reco::ClusterParameters& clusterParams1 = outputClusterList.back();
693 
694  reco::PrincipalComponents& fullPCA(clusterToBreak.getFullPCA());
695  Eigen::Vector3f fullPrimaryVec(fullPCA.getEigenVectors().row(2));
696 
697  if (makeCandidateCluster(fullPrimaryVec, clusterParams1, hitList.begin(), kinkItr, level))
698  {
699  outputClusterList.push_back(reco::ClusterParameters());
700 
701  reco::ClusterParameters& clusterParams2 = outputClusterList.back();
702 
703  makeCandidateCluster(fullPrimaryVec, clusterParams2, kinkItr, hitList.end(), level);
704 
705  if (fFillHistograms)
706  {
707  fillConvexHullHists(clusterParams1, false);
708  fillConvexHullHists(clusterParams2, false);
709  }
710  }
711 
712  // If we did not make 2 clusters then be sure to clear the output list
713  if (outputClusterList.size() != 2) outputClusterList.clear();
714  }
715 
716  return !outputClusterList.empty();
717 }
718 
720 {
721  // Set up container to keep track of edges
722  KinkTupleVec kinkTupleVec;
723 
724  reco::ConvexHull& convexHull = clusterToBreak.getConvexHull();
725  reco::ProjectedPointList& pointList = convexHull.getProjectedPointList();
726  reco::ConvexHullKinkTupleList& kinkPointList = convexHull.getConvexHullKinkPoints();
727 
728  for(auto& kink : kinkPointList)
729  {
730  // Make an instance of the vec value to avoid copying if we can...
731  kinkTupleVec.push_back(KinkTuple());
732 
733  KinkTuple& kinkTuple = kinkTupleVec.back();
734 
735  std::get<1>(kinkTuple) = kink;
736 
737  // Recover vectors, want them pointing away from intersection point
738  Eigen::Vector2f firstEdge = -std::get<1>(kink);
739  HitOrderTupleList& firstList = std::get<2>(kinkTuple);
740  HitOrderTupleList& secondList = std::get<3>(kinkTuple);
741 
742  orderHitsAlongEdge(pointList, std::get<0>(kink), firstEdge, firstList);
743 
744  if (firstList.size() > fMinTinyClusterSize)
745  {
746  Eigen::Vector2f secondEdge = std::get<2>(kink);
747 
748  orderHitsAlongEdge(pointList, std::get<0>(kink), secondEdge, secondList);
749 
750  if (secondList.size() > fMinTinyClusterSize)
751  std::get<0>(kinkTuple) = std::min(firstList.size(),secondList.size());
752  }
753 
754  // Special handling...
755  if (firstList.size() + secondList.size() > pointList.size())
756  {
757  if (firstList.size() > secondList.size()) pruneHitOrderTupleLists(firstList,secondList);
758  else pruneHitOrderTupleLists(secondList,firstList);
759 
760  std::get<0>(kinkTuple) = std::min(firstList.size(),secondList.size());
761  }
762 
763  if (std::get<0>(kinkTuple) < int(fMinTinyClusterSize)) kinkTupleVec.pop_back();
764  }
765 
766  // No work if the list is empty
767  if (!kinkTupleVec.empty())
768  {
769  // If more than one then want the kink with the most elements both sizes
770  std::sort(kinkTupleVec.begin(),kinkTupleVec.end(),[](const auto& left,const auto& right){return std::get<0>(left) > std::get<0>(right);});
771 
772  // Recover the kink point
773  KinkTuple& kinkTuple = kinkTupleVec.front();
774 
775  // Set up to split the input cluster
776  outputClusterList.push_back(reco::ClusterParameters());
777 
778  reco::ClusterParameters& clusterParams1 = outputClusterList.back();
779 
780  reco::PrincipalComponents& fullPCA(clusterToBreak.getFullPCA());
781  Eigen::Vector3f fullPrimaryVec(fullPCA.getEigenVectors().row(2));
782 
783  if (makeCandidateCluster(fullPrimaryVec, clusterParams1, std::get<2>(kinkTuple), level))
784  {
785  outputClusterList.push_back(reco::ClusterParameters());
786 
787  reco::ClusterParameters& clusterParams2 = outputClusterList.back();
788 
789  makeCandidateCluster(fullPrimaryVec, clusterParams2, std::get<3>(kinkTuple), level);
790 
791  if (fFillHistograms)
792  {
793  fillConvexHullHists(clusterParams1, false);
794  fillConvexHullHists(clusterParams2, false);
795  }
796  }
797 
798  // If we did not make 2 clusters then be sure to clear the output list
799  if (outputClusterList.size() != 2) outputClusterList.clear();
800  }
801 
802  return !outputClusterList.empty();
803 }
804 
806  const reco::ProjectedPoint& point,
807  const Eigen::Vector2f& edge,
808  HitOrderTupleList& orderedList) const
809 {
810  // Use the input kink point as the start point of the edge
811  Eigen::Vector2f kinkPos(std::get<0>(point),std::get<1>(point));
812 
813  // Loop over the input hits
814  for (const auto& hit : hitList)
815  {
816  // Now we need to calculate the doca and poca...
817  // Start by getting this hits position
818  Eigen::Vector2f hitPos(std::get<0>(hit),std::get<1>(hit));
819 
820  // Form a TVector from this to the cluster average position
821  Eigen::Vector2f hitToKinkVec = hitPos - kinkPos;
822 
823  // With this we can get the arclength to the doca point
824  float arcLenToPoca = hitToKinkVec.dot(edge);
825 
826  // Require the hit to not be past the kink point
827  if (arcLenToPoca < 0.) continue;
828 
829  // Get the coordinates along the axis for this point
830  Eigen::Vector2f pocaPos = kinkPos + arcLenToPoca * edge;
831 
832  // Now get doca and poca
833  Eigen::Vector2f pocaPosToHitPos = hitPos - pocaPos;
834  float pocaToAxis = pocaPosToHitPos.norm();
835 
836  std::cout << "-- arcLenToPoca: " << arcLenToPoca << ", doca: " << pocaToAxis << std::endl;
837 
838  orderedList.emplace_back(arcLenToPoca,pocaToAxis,hit);
839  }
840 
841  // Sort the list in order of ascending distance from the kink point
842  orderedList.sort([](const auto& left,const auto& right){return std::get<0>(left) < std::get<0>(right);});
843 
844  return;
845 }
846 
848 {
849  // Assume the first list is the short one, so we loop through the elements of that list..
850  HitOrderTupleList::iterator shortItr = shortList.begin();
851 
852  while(shortItr != shortList.end())
853  {
854  // Recover the search key
855  const reco::ClusterHit3D* hit3D = std::get<2>(std::get<2>(*shortItr));
856 
857  // Ok, find the corresponding point in the other list...
858  HitOrderTupleList::iterator longItr = std::find_if(longList.begin(),longList.end(),[&hit3D](const auto& elem){return hit3D == std::get<2>(std::get<2>(elem));});
859 
860  if (longItr != longList.end())
861  {
862  if (std::get<1>(*longItr) < std::get<1>(*shortItr))
863  {
864  shortItr = shortList.erase(shortItr);
865  }
866  else
867  {
868  longItr = longList.erase(longItr);
869  shortItr++;
870  }
871  }
872  else shortItr++;
873  }
874 
875 
876  return;
877 }
878 
880 {
881  // Set up container to keep track of edges
882  using DistEdgeTuple = std::tuple<float, const reco::EdgeTuple*>;
883  using DistEdgeTupleVec = std::vector<DistEdgeTuple>;
884 
885  DistEdgeTupleVec distEdgeTupleVec;
886 
887  reco::ProjectedPointList::const_iterator extremePointListItr = clusterToBreak.getConvexHull().getConvexHullExtremePoints().begin();
888 
889  const reco::ClusterHit3D* firstEdgeHit = std::get<2>(*extremePointListItr++);
890  const reco::ClusterHit3D* secondEdgeHit = std::get<2>(*extremePointListItr);
891  Eigen::Vector3f edgeVec(secondEdgeHit->getPosition()[0] - firstEdgeHit->getPosition()[0],
892  secondEdgeHit->getPosition()[1] - firstEdgeHit->getPosition()[1],
893  secondEdgeHit->getPosition()[2] - firstEdgeHit->getPosition()[2]);
894  double edgeLen = edgeVec.norm();
895 
896  // normalize it
897  edgeVec.normalize();
898 
899  // Now loop through all the edges and search for the furthers point
900  for(const auto& edge : clusterToBreak.getBestEdgeList())
901  {
902  const reco::ClusterHit3D* nextEdgeHit = std::get<0>(edge); // recover the first point
903 
904  // Create vector to this point from the longest edge
905  Eigen::Vector3f hitToEdgeVec(nextEdgeHit->getPosition()[0] - firstEdgeHit->getPosition()[0],
906  nextEdgeHit->getPosition()[1] - firstEdgeHit->getPosition()[1],
907  nextEdgeHit->getPosition()[2] - firstEdgeHit->getPosition()[2]);
908 
909  // Get projection
910  float hitProjection = hitToEdgeVec.dot(edgeVec);
911 
912  // Require that the point is really "opposite" the longest edge
913  if (hitProjection > 0. && hitProjection < edgeLen)
914  {
915  Eigen::Vector3f distToHitVec = hitToEdgeVec - hitProjection * edgeVec;
916  float distToHit = distToHitVec.norm();
917 
918  distEdgeTupleVec.emplace_back(distToHit,&edge);
919  }
920  }
921 
922  std::sort(distEdgeTupleVec.begin(),distEdgeTupleVec.end(),[](const auto& left,const auto& right){return std::get<0>(left) > std::get<0>(right);});
923 
924  reco::PrincipalComponents& fullPCA(clusterToBreak.getFullPCA());
925  Eigen::Vector3f fullPrimaryVec(fullPCA.getEigenVectors().row(2));
926 
927  // Recover our hits
928  reco::HitPairListPtr& hitList = clusterToBreak.getHitPairListPtr();
929 
930  // Calculate the doca to the PCA primary axis for each 3D hit
931  // Importantly, this also gives us the arclength along that axis to the hit
932  fPCAAlg.PCAAnalysis_calc3DDocas(hitList, fullPCA);
933 
934  // Sort the hits along the PCA
935  hitList.sort([](const auto& left, const auto& right){return left->getArclenToPoca() < right->getArclenToPoca();});
936 
937  // Get a temporary container to hol
938  float usedDefectDist(0.);
939 
940  for(const auto& distEdgeTuple : distEdgeTupleVec)
941  {
942  const reco::EdgeTuple& edgeTuple = *std::get<1>(distEdgeTuple);
943  const reco::ClusterHit3D* edgeHit = std::get<0>(edgeTuple);
944 
945  usedDefectDist = std::get<0>(distEdgeTuple);
946 
947  // Now find the hit identified above as furthest away
948  reco::HitPairListPtr::iterator vertexItr = std::find(hitList.begin(),hitList.end(),edgeHit);
949 
950  // Make sure enough hits either side, otherwise we just keep the current cluster
951  if (vertexItr == hitList.end() || std::distance(hitList.begin(),vertexItr) < int(fMinTinyClusterSize) || std::distance(vertexItr,hitList.end()) < int(fMinTinyClusterSize)) continue;
952 
953  outputClusterList.push_back(reco::ClusterParameters());
954 
955  reco::ClusterParameters& clusterParams1 = outputClusterList.back();
956 
957  if (makeCandidateCluster(fullPrimaryVec, clusterParams1, hitList.begin(), vertexItr, level))
958  {
959  outputClusterList.push_back(reco::ClusterParameters());
960 
961  reco::ClusterParameters& clusterParams2 = outputClusterList.back();
962 
963  if (makeCandidateCluster(fullPrimaryVec, clusterParams2, vertexItr, hitList.end(), level))
964  {
965  if (fFillHistograms)
966  {
967  fSubMaxDefect->Fill(std::get<0>(distEdgeTupleVec.front()), 1.);
968  fSubUsedDefect->Fill(usedDefectDist, 1.);
969  }
970  break;
971  }
972  }
973 
974  // If here then we could not make two valid clusters and so we try again
975  outputClusterList.clear();
976  }
977 
978  return !outputClusterList.empty();
979 }
980 
981 bool ConvexHullPathFinder::breakClusterInHalf(reco::ClusterParameters& clusterToBreak, reco::ClusterParametersList& outputClusterList, int level) const
982 {
983  reco::PrincipalComponents& fullPCA(clusterToBreak.getFullPCA());
984  Eigen::Vector3f fullPrimaryVec(fullPCA.getEigenVectors().row(2));
985 
986  // Recover our hits
987  reco::HitPairListPtr& hitList = clusterToBreak.getHitPairListPtr();
988 
989  // Calculate the doca to the PCA primary axis for each 3D hit
990  // Importantly, this also gives us the arclength along that axis to the hit
991  fPCAAlg.PCAAnalysis_calc3DDocas(hitList, fullPCA);
992 
993  // Sort the hits along the PCA
994  hitList.sort([](const auto& left, const auto& right){return left->getArclenToPoca() < right->getArclenToPoca();});
995 
996  reco::HitPairListPtr::iterator vertexItr = hitList.begin();
997 
998  std::advance(vertexItr, hitList.size()/2);
999 
1000  outputClusterList.push_back(reco::ClusterParameters());
1001 
1002  reco::ClusterParameters& clusterParams1 = outputClusterList.back();
1003 
1004  if (makeCandidateCluster(fullPrimaryVec, clusterParams1, hitList.begin(), vertexItr, level))
1005  {
1006  outputClusterList.push_back(reco::ClusterParameters());
1007 
1008  reco::ClusterParameters& clusterParams2 = outputClusterList.back();
1009 
1010  makeCandidateCluster(fullPrimaryVec, clusterParams2, vertexItr, hitList.end(), level);
1011  }
1012 
1013  if (outputClusterList.size() != 2) outputClusterList.clear();
1014 
1015  return !outputClusterList.empty();
1016 }
1017 
1019 {
1020  // Idea here is to scan the input hit list (assumed ordered along the current PCA) and look for "large" gaps
1021  // Here a gap is determined when the hits were ordered by their distance along the primary PCA to their doca to it.
1022  reco::PrincipalComponents& fullPCA(clusterToBreak.getFullPCA());
1023 
1024  // Recover our hits
1025  reco::HitPairListPtr& hitList = clusterToBreak.getHitPairListPtr();
1026 
1027  // Calculate the doca to the PCA primary axis for each 3D hit
1028  // Importantly, this also gives us the arclength along that axis to the hit
1029  fPCAAlg.PCAAnalysis_calc3DDocas(hitList, fullPCA);
1030 
1031  // Sort the hits along the PCA
1032  hitList.sort([](const auto& left, const auto& right){return left->getArclenToPoca() < right->getArclenToPoca();});
1033 
1034  // Loop through the input hit list and keep track of first hit of largest gap
1035  reco::HitPairListPtr::iterator bigGapHitItr = hitList.begin();
1036  float biggestGap = 0.;
1037 
1038  reco::HitPairListPtr::iterator lastHitItr = hitList.begin();
1039 
1040  for(reco::HitPairListPtr::iterator hitItr = hitList.begin(); hitItr != hitList.end(); hitItr++)
1041  {
1042  float currentGap = std::abs((*hitItr)->getArclenToPoca() - (*lastHitItr)->getArclenToPoca());
1043 
1044  if (currentGap > biggestGap)
1045  {
1046  bigGapHitItr = hitItr; // Note that this is an iterator and will be the "end" going from begin, and "begin" for second half
1047  biggestGap = currentGap;
1048  }
1049 
1050  lastHitItr = hitItr;
1051  }
1052 
1053  // Require some minimum gap size...
1054  if (biggestGap > fMinGapSize)
1055  {
1056  outputClusterList.push_back(reco::ClusterParameters());
1057 
1058  reco::ClusterParameters& clusterParams1 = outputClusterList.back();
1059 
1060  reco::PrincipalComponents& fullPCA(clusterToBreak.getFullPCA());
1061  Eigen::Vector3f fullPrimaryVec(fullPCA.getEigenVectors().row(2));
1062 
1063  if (makeCandidateCluster(fullPrimaryVec, clusterParams1, hitList.begin(), bigGapHitItr, level))
1064  {
1065  outputClusterList.push_back(reco::ClusterParameters());
1066 
1067  reco::ClusterParameters& clusterParams2 = outputClusterList.back();
1068 
1069  makeCandidateCluster(fullPrimaryVec, clusterParams2, bigGapHitItr, hitList.end(), level);
1070  }
1071 
1072  if (outputClusterList.size() != 2) outputClusterList.clear();
1073  }
1074 
1075  return !outputClusterList.empty();
1076 }
1077 
1078 
1079 void ConvexHullPathFinder::buildConvexHull(reco::ClusterParameters& clusterParameters, int level) const
1080 {
1081  // set an indention string
1082  std::string minuses(level/2, '-');
1083  std::string indent(level/2, ' ');
1084 
1085  indent += minuses;
1086 
1087  // The plan is to build the enclosing 2D polygon around the points in the PCA plane of most spread for this cluster
1088  // To do so we need to start by building a list of 2D projections onto the plane of most spread...
1089  reco::PrincipalComponents& pca = clusterParameters.getFullPCA();
1090 
1091  // Recover the parameters from the Principal Components Analysis that we need to project and accumulate
1092  const Eigen::Vector3f& pcaCenter = pca.getAvePosition();
1093  reco::ConvexHull& convexHull = clusterParameters.getConvexHull();
1094  reco::ProjectedPointList& pointList = convexHull.getProjectedPointList();
1095 
1096  // Loop through hits and do projection to plane
1097  for(const auto& hit3D : clusterParameters.getHitPairListPtr())
1098  {
1099  Eigen::Vector3f pcaToHitVec(hit3D->getPosition()[0] - pcaCenter(0),
1100  hit3D->getPosition()[1] - pcaCenter(1),
1101  hit3D->getPosition()[2] - pcaCenter(2));
1102  Eigen::Vector3f pcaToHit = pca.getEigenVectors() * pcaToHitVec;
1103 
1104  pointList.emplace_back(pcaToHit(1),pcaToHit(2),hit3D);
1105  }
1106 
1107  // Sort the point vec by increasing x, then increase y
1108  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);});
1109 
1110  // containers for finding the "best" hull...
1111  std::vector<ConvexHull> convexHullVec;
1112  std::vector<reco::ProjectedPointList> rejectedListVec;
1113  bool increaseDepth(pointList.size() > 3);
1114  float lastArea(std::numeric_limits<float>::max());
1115 
1116  while(increaseDepth)
1117  {
1118  // Get another convexHull container
1119  convexHullVec.push_back(ConvexHull(pointList, fConvexHullKinkAngle, fConvexHullMinSep));
1120  rejectedListVec.push_back(reco::ProjectedPointList());
1121 
1122  const ConvexHull& convexHull = convexHullVec.back();
1123  reco::ProjectedPointList& rejectedList = rejectedListVec.back();
1124  const reco::ProjectedPointList& convexHullPoints = convexHull.getConvexHull();
1125 
1126  increaseDepth = false;
1127 
1128  if (convexHull.getConvexHullArea() > 0.)
1129  {
1130  if (convexHullVec.size() < 2 || convexHull.getConvexHullArea() < 0.8 * lastArea)
1131  {
1132  for(auto& point : convexHullPoints)
1133  {
1134  pointList.remove(point);
1135  rejectedList.emplace_back(point);
1136  }
1137  lastArea = convexHull.getConvexHullArea();
1138 // increaseDepth = true;
1139  }
1140  }
1141  }
1142 
1143  // do we have a valid convexHull?
1144  while(!convexHullVec.empty() && convexHullVec.back().getConvexHullArea() < 0.5)
1145  {
1146  convexHullVec.pop_back();
1147  rejectedListVec.pop_back();
1148  }
1149 
1150  // If we found the convex hull then build edges around the region
1151  if (!convexHullVec.empty())
1152  {
1153  size_t nRejectedTotal(0);
1154  reco::HitPairListPtr hitPairListPtr = clusterParameters.getHitPairListPtr();
1155 
1156  for(const auto& rejectedList : rejectedListVec)
1157  {
1158  nRejectedTotal += rejectedList.size();
1159 
1160  for(const auto& rejectedPoint : rejectedList)
1161  {
1162  if (convexHullVec.back().findNearestDistance(rejectedPoint) > 0.5)
1163  hitPairListPtr.remove(std::get<2>(rejectedPoint));
1164  }
1165  }
1166 
1167  // Now add "edges" to the cluster to describe the convex hull (for the display)
1168  reco::ProjectedPointList& convexHullPointList = convexHull.getConvexHullPointList();
1169  reco::Hit3DToEdgeMap& edgeMap = convexHull.getConvexHullEdgeMap();
1170  reco::EdgeList& edgeList = convexHull.getConvexHullEdgeList();
1171 
1172  reco::ProjectedPoint lastPoint = convexHullVec.back().getConvexHull().front();
1173 
1174  for(auto& curPoint : convexHullVec.back().getConvexHull())
1175  {
1176  if (curPoint == lastPoint) continue;
1177 
1178  const reco::ClusterHit3D* lastPoint3D = std::get<2>(lastPoint);
1179  const reco::ClusterHit3D* curPoint3D = std::get<2>(curPoint);
1180 
1181  float distBetweenPoints = (curPoint3D->getPosition()[0] - lastPoint3D->getPosition()[0]) * (curPoint3D->getPosition()[0] - lastPoint3D->getPosition()[0])
1182  + (curPoint3D->getPosition()[1] - lastPoint3D->getPosition()[1]) * (curPoint3D->getPosition()[1] - lastPoint3D->getPosition()[1])
1183  + (curPoint3D->getPosition()[2] - lastPoint3D->getPosition()[2]) * (curPoint3D->getPosition()[2] - lastPoint3D->getPosition()[2]);
1184 
1185  distBetweenPoints = std::sqrt(distBetweenPoints);
1186 
1187  reco::EdgeTuple edge(lastPoint3D,curPoint3D,distBetweenPoints);
1188 
1189  convexHullPointList.push_back(curPoint);
1190  edgeMap[lastPoint3D].push_back(edge);
1191  edgeMap[curPoint3D].push_back(edge);
1192  edgeList.emplace_back(edge);
1193 
1194  lastPoint = curPoint;
1195  }
1196 
1197  // Store the "extreme" points
1198  const ConvexHull::PointList& extremePoints = convexHullVec.back().getExtremePoints();
1199  reco::ProjectedPointList& extremePointList = convexHull.getConvexHullExtremePoints();
1200 
1201  for(const auto& point : extremePoints) extremePointList.push_back(point);
1202 
1203  // Store the "kink" points
1204  const reco::ConvexHullKinkTupleList& kinkPoints = convexHullVec.back().getKinkPoints();
1205  reco::ConvexHullKinkTupleList& kinkPointList = convexHull.getConvexHullKinkPoints();
1206 
1207  for(const auto& kink : kinkPoints) kinkPointList.push_back(kink);
1208  }
1209 
1210  return;
1211 }
1212 
1214 {
1215  reco::ProjectedPointList& convexHullPoints = clusterParameters.getConvexHull().getConvexHullPointList();
1216 
1217  if (convexHullPoints.size() > 2)
1218  {
1219  reco::ProjectedPointList::iterator pointItr = convexHullPoints.begin();
1220 
1221  // Advance to the second to last element
1222  std::advance(pointItr, convexHullPoints.size() - 2);
1223 
1224  reco::ProjectedPoint lastPoint = *pointItr++;
1225 
1226  // Reset pointer to the first element
1227  pointItr = convexHullPoints.begin();
1228 
1229  reco::ProjectedPoint curPoint = *pointItr++;
1230  Eigen::Vector2f lastEdge(std::get<0>(curPoint) - std::get<0>(lastPoint), std::get<1>(curPoint) - std::get<1>(lastPoint));
1231 
1232  lastEdge.normalize();
1233 
1234  while(pointItr != convexHullPoints.end())
1235  {
1236  reco::ProjectedPoint& nextPoint = *pointItr++;
1237 
1238  Eigen::Vector2f nextEdge(std::get<0>(nextPoint) - std::get<0>(curPoint), std::get<1>(nextPoint) - std::get<1>(curPoint));
1239  float nextEdgeLen = nextEdge.norm();
1240 
1241  nextEdge.normalize();
1242 
1243  float cosLastNextEdge = lastEdge.dot(nextEdge);
1244 
1245  if (top)
1246  {
1247  fTopConvexCosEdge->Fill(cosLastNextEdge, 1.);
1248  fTopConvexEdgeLen->Fill(std::min(nextEdgeLen,float(49.9)), 1.);
1249  }
1250  else
1251  {
1252  fSubConvexCosEdge->Fill(cosLastNextEdge, 1.);
1253  fSubConvexEdgeLen->Fill(std::min(nextEdgeLen,float(49.9)), 1.);
1254  }
1255 
1256  if (nextEdgeLen > fConvexHullMinSep) lastEdge = nextEdge;
1257 
1258  curPoint = nextPoint;
1259  }
1260  }
1261 
1262  return;
1263 }
1264 
1265 float ConvexHullPathFinder::closestApproach(const Eigen::Vector3f& P0,
1266  const Eigen::Vector3f& u0,
1267  const Eigen::Vector3f& P1,
1268  const Eigen::Vector3f& u1,
1269  Eigen::Vector3f& poca0,
1270  Eigen::Vector3f& poca1) const
1271 {
1272  // Technique is to compute the arclength to each point of closest approach
1273  Eigen::Vector3f w0 = P0 - P1;
1274  float a(1.);
1275  float b(u0.dot(u1));
1276  float c(1.);
1277  float d(u0.dot(w0));
1278  float e(u1.dot(w0));
1279  float den(a * c - b * b);
1280 
1281  float arcLen0 = (b * e - c * d) / den;
1282  float arcLen1 = (a * e - b * d) / den;
1283 
1284  poca0 = P0 + arcLen0 * u0;
1285  poca1 = P1 + arcLen1 * u1;
1286 
1287  return (poca0 - poca1).norm();
1288 }
1289 
1290 float ConvexHullPathFinder::findConvexHullEndPoints(const reco::EdgeList& convexHull, const reco::ClusterHit3D* first3D, const reco::ClusterHit3D* last3D) const
1291 {
1292  float largestDistance(0.);
1293 
1294  // Note that edges are vectors and that the convex hull edge list will be ordered
1295  // 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
1296  // meaning that the dot product of the two edges becomes negative.
1297  reco::EdgeList::const_iterator firstEdgeItr = convexHull.begin();
1298 
1299  while(firstEdgeItr != convexHull.end())
1300  {
1301  reco::EdgeList::const_iterator nextEdgeItr = firstEdgeItr;
1302 
1303 // Eigen::Vector2f firstEdgeVec(std::get<3>(*firstEdgeItr),std::get<);
1304 // Eigen::Vector2f lastPrimaryVec(lastPCA.getEigenVectors()[0][0],lastPCA.getEigenVectors()[0][1],lastPCA.getEigenVectors()[0][2]);
1305 // float cosToLast = newPrimaryVec.dot(lastPrimaryVec);
1306 
1307  while(++nextEdgeItr != convexHull.end())
1308  {
1309 
1310  }
1311  }
1312 
1313  return largestDistance;
1314 }
1315 
1316 DEFINE_ART_CLASS_TOOL(ConvexHullPathFinder)
1317 } // namespace lar_cluster3d
size_t fMinTinyClusterSize
Minimum size for a &quot;tiny&quot; cluster.
float fConvexHullMinSep
Min hit separation to conisder in convex hull.
void PCAAnalysis_calc3DDocas(const reco::HitPairListPtr &hitPairVector, const reco::PrincipalComponents &pca) const
bool getSvdOK() const
Definition: Cluster3D.h:243
std::tuple< int, reco::ConvexHullKinkTuple, HitOrderTupleList, HitOrderTupleList > KinkTuple
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
bool breakClusterByMaxDefect(reco::ClusterParameters &, reco::ClusterParametersList &, int level) const
const Eigen::Vector3f getPosition() const
Definition: Cluster3D.h:157
std::tuple< float, float, reco::ProjectedPoint > HitOrderTuple
std::list< ProjectedPoint > ProjectedPointList
Definition: Cluster3D.h:352
walls no right
Definition: selectors.fcl:105
reco::PrincipalComponents & getSkeletonPCA()
Definition: Cluster3D.h:477
std::list< Point > PointList
The list of the projected points.
Definition: ConvexHull.h:34
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.
reco::EdgeList & getBestEdgeList()
Definition: Cluster3D.h:480
reco::HitPairListPtr & getHitPairListPtr()
Definition: Cluster3D.h:475
process_name hit
Definition: cheaterreco.fcl:51
void configure(fhicl::ParameterSet const &pset) override
Interface for configuring the particular algorithm tool.
IClusterModAlg interface class definiton.
process_name gaushit a
walls no top
Definition: selectors.fcl:105
void ModifyClusters(reco::ClusterParametersList &) const override
Scan an input collection of clusters and modify those according to the specific implementing algorith...
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
double distance(geo::Point_t const &point, CathodeDesc_t const &cathode)
Returns the distance of a point from the cathode.
void fillConvexHullHists(reco::ClusterParameters &, bool) const
reco::PrincipalComponents & getFullPCA()
Definition: Cluster3D.h:476
reco::Hit3DToEdgeMap & getConvexHullEdgeMap()
Definition: Cluster3D.h:384
reco::ProjectedPointList & getConvexHullExtremePoints()
Definition: Cluster3D.h:386
reco::ConvexHullKinkTupleList & getConvexHullKinkPoints()
Definition: Cluster3D.h:387
Define a container for working with the convex hull.
Definition: Cluster3D.h:359
const Eigen::Vector3f & getAvePosition() const
Definition: Cluster3D.h:247
void pruneHitOrderTupleLists(HitOrderTupleList &, HitOrderTupleList &) const
ConvexHullPathFinder(const fhicl::ParameterSet &)
Constructor.
std::list< ConvexHullKinkTuple > ConvexHullKinkTupleList
Definition: Cluster3D.h:354
std::tuple< const reco::ClusterHit3D *, const reco::ClusterHit3D *, double > EdgeTuple
Definition: Cluster3D.h:343
std::list< const reco::ClusterHit3D * > HitPairListPtr
Definition: Cluster3D.h:334
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 ...
void buildConvexHull(reco::ClusterParameters &clusterParameters, int level=0) const
std::unique_ptr< lar_cluster3d::IClusterAlg > fClusterAlg
Tools.
bool breakClusterByKinksTrial(reco::ClusterParameters &, reco::ClusterParametersList &, int level) const
bool makeCandidateCluster(Eigen::Vector3f &, reco::ClusterParameters &, reco::HitPairListPtr::iterator, reco::HitPairListPtr::iterator, int) const
reco::ConvexHull & getConvexHull()
Definition: Cluster3D.h:481
void initializeHistograms(art::TFileDirectory &) override
Interface for initializing histograms if they are desired Note that the idea is to put hisgtograms in...
bool breakClusterByKinks(reco::ClusterParameters &, reco::ClusterParametersList &, int level) const
float fMinEigen0To1Ratio
Minimum ratio of eigen 0 to 1 to continue breaking.
float closestApproach(const Eigen::Vector3f &, const Eigen::Vector3f &, const Eigen::Vector3f &, const Eigen::Vector3f &, Eigen::Vector3f &, Eigen::Vector3f &) const
void orderHitsAlongEdge(const reco::ProjectedPointList &, const reco::ProjectedPoint &, const Eigen::Vector2f &, HitOrderTupleList &) const
float fConvexHullKinkAngle
Angle to declare a kink in convex hull calc.
This provides an art tool interface definition for 3D Cluster algorithms.
bool breakClusterAtBigGap(reco::ClusterParameters &, reco::ClusterParametersList &, int level) const
do i e
T copy(T const &v)
bool breakClusterInHalf(reco::ClusterParameters &, reco::ClusterParametersList &, int level) const
std::unordered_map< const reco::ClusterHit3D *, reco::EdgeList > Hit3DToEdgeMap
Definition: Cluster3D.h:346
reco::ProjectedPointList & getConvexHullPointList()
Definition: Cluster3D.h:383
process_name physics producers generator hPHist_pi physics producers generator P0
reco::ProjectedPointList & getProjectedPointList()
Definition: Cluster3D.h:382
std::pair< MinMaxPoints, MinMaxPoints > MinMaxPointPair
reco::EdgeList & getConvexHullEdgeList()
Definition: Cluster3D.h:385
float findConvexHullEndPoints(const reco::EdgeList &, const reco::ClusterHit3D *, const reco::ClusterHit3D *) const
std::tuple< float, float, const reco::ClusterHit3D * > ProjectedPoint
Projected coordinates and pointer to hit.
Definition: Cluster3D.h:351
BEGIN_PROLOG don t mess with this pandoraTrackGausCryoW true
float fMinGapSize
Minimum gap size to break at gaps.
BEGIN_PROLOG could also be cout
float getTimeToExecute() const override
If monitoring, recover the time to execute a particular function.
std::list< ClusterParameters > ClusterParametersList
Definition: Cluster3D.h:403
const EigenVectors & getEigenVectors() const
Definition: Cluster3D.h:246
std::pair< reco::ProjectedPoint, reco::ProjectedPoint > MinMaxPoints
bool completeCandidateCluster(Eigen::Vector3f &, reco::ClusterParameters &, int) const