HoughSeedFinderAlg.cxx
Go to the documentation of this file.
1 /**
2  * @file HoughSeedFinderAlg.cxx
3  *
4  * @brief Implementation of the Seed Finder Algorithm based on a Hough Transform
5  *
6  * The intent of this algorithm is to take an input list of 3D space points and from those
7  * to find candidate track start points and directions. It does so by first performing a
8  * Principal Components Analysis of the input 3D hits and then projects them to the plane
9  * of largest spread. A standard Hough Transform method is then applied to attempt to identify
10  * straight line segments which can be used as seeds to the kalman filter tracker.
11  */
12 
13 // The main include
15 // Framework Includes
16 
17 // LArSoft includes
20 
21 // ROOT includes
22 #include "TCanvas.h"
23 #include "TDecompSVD.h"
24 #include "TFrame.h"
25 #include "TH2D.h"
26 #include "TMatrixD.h"
27 #include "TVectorD.h"
28 
29 // std includes
30 #include <memory>
31 
32 //------------------------------------------------------------------------------------------------------------------------------------------
33 
34 //------------------------------------------------------------------------------------------------------------------------------------------
35 // implementation follows
36 
37 namespace lar_cluster3d {
38 
40  : m_minimum3DHits(5)
41  , m_thetaBins(360)
42  , m_rhoBins(75)
43  , m_hiThresholdMin(5)
44  , m_hiThresholdFrac(.05)
45  , m_loThresholdFrac(0.85)
46  , m_numSeed2DHits(80)
47  , m_numAveDocas(6.)
48  , m_numSkippedHits(10)
49  , m_maxLoopsPerCluster(3)
50  , m_maximumGap(5.)
51  , m_pcaAlg(pset.get<fhicl::ParameterSet>("PrincipalComponentsAlg"))
52  , m_displayHist(false)
53  {
54  m_minimum3DHits = pset.get<size_t>("Minimum3DHits", 5);
55  m_thetaBins = pset.get<int>("ThetaBins", 360);
56  m_rhoBins = pset.get<int>("HalfRhoBins", 75);
57  m_hiThresholdMin = pset.get<size_t>("HiThresholdMin", 5);
58  m_hiThresholdFrac = pset.get<double>("HiThresholdFrac", 0.05);
59  m_loThresholdFrac = pset.get<double>("LoThresholdFrac", 0.85);
60  m_numSeed2DHits = pset.get<size_t>("NumSeed2DHits", 80);
61  m_numAveDocas = pset.get<double>("NumAveDocas", 6.);
62  m_numSkippedHits = pset.get<int>("NumSkippedHits", 10);
63  m_maxLoopsPerCluster = pset.get<int>("MaxLoopsPerCluster", 3);
64  m_maximumGap = pset.get<double>("MaximumGap", 5.);
65  m_displayHist = pset.get<bool>("DisplayHoughHist", false);
67  }
68 
69  //------------------------------------------------------------------------------------------------------------------------------------------
70 
71  class AccumulatorValues {
72  /**
73  * @brief A utility class to contain the values of a given "bin" in Hough Space
74  *
75  * Specifically, this is keeping track of the projected x,y coordinates of a given
76  * 3D hit projected to the plane of largest spread in PCA and an interator to
77  * that hit in the input container
78  */
79  public:
80  AccumulatorValues() : m_position(0., 0., 0.) {}
81  AccumulatorValues(const Eigen::Vector3f& position,
83  : m_position(position), m_hit3DIterator(itr)
84  {}
85 
86  const Eigen::Vector3f&
87  getPosition() const
88  {
89  return m_position;
90  }
93  {
94  return m_hit3DIterator;
95  }
96 
97  private:
98  Eigen::Vector3f
99  m_position; ///< We really only need the x,y coordinates here but keep all three for now
101  m_hit3DIterator; ///< This will be used to take us back to our 3D hit
102  };
103 
104  typedef std::vector<AccumulatorValues> AccumulatorValuesVec;
105 
107  /**
108  * @brief A utility class used to accumulate the above values
109  *
110  * One of these objects will exist for each "bin" in rho-theta space and this will
111  * be used to accumulate the 3D hits which contribute to this bin
112  */
113  public:
114  AccumulatorBin() : m_visited(false), m_noise(false), m_inCluster(false) {}
115  AccumulatorBin(AccumulatorValues& values) : m_visited(false), m_noise(false), m_inCluster(false)
116  {
117  m_accumulatorValuesVec.push_back(values);
118  }
119 
120  void
122  {
123  m_visited = true;
124  }
125  void
127  {
128  m_noise = true;
129  }
130  void
132  {
133  m_inCluster = true;
134  }
135 
136  void
138  {
139  m_accumulatorValuesVec.push_back(value);
140  }
141 
142  bool
143  isVisited() const
144  {
145  return m_visited;
146  }
147  bool
148  isNoise() const
149  {
150  return m_noise;
151  }
152  bool
153  isInCluster() const
154  {
155  return m_inCluster;
156  }
157 
158  const AccumulatorValuesVec&
160  {
161  return m_accumulatorValuesVec;
162  }
163 
164  private:
165  bool m_visited;
166  bool m_noise;
167  bool m_inCluster;
168  AccumulatorValuesVec m_accumulatorValuesVec;
169  };
170 
172  /**
173  * @brief This is used to sort "Hough Clusters" by the maximum entries in a bin
174  */
175  public:
177  {}
178 
179  bool
182  {
183  size_t peakCountLeft(0);
184  size_t peakCountRight(0);
185 
186  for (const auto& binIndex : left)
187  peakCountLeft = std::max(peakCountLeft, m_accMap[binIndex].getAccumulatorValues().size());
188  for (const auto& binIndex : right)
189  peakCountRight = std::max(peakCountRight, m_accMap[binIndex].getAccumulatorValues().size());
190 
191  return peakCountLeft > peakCountRight;
192  }
193 
194  private:
196  };
197 
199  /**
200  * @brief This is used to sort bins in Hough Space
201  */
202  bool
205  {
206  size_t leftSize = left->second.getAccumulatorValues().size();
207  size_t rightSize = right->second.getAccumulatorValues().size();
208 
209  return leftSize > rightSize;
210  }
211  };
212 
213  void
215  RhoThetaAccumulatorBinMap& rhoThetaAccumulatorBinMap,
216  HoughCluster& neighborPts,
217  size_t threshold) const
218  {
219  /**
220  * @brief Does a query of nearest neighbors to look for matching bins
221  */
222 
223  // We simply loop over the nearest indices and see if we have any friends over threshold
224  for (int rhoIdx = curBin.first - 1; rhoIdx <= curBin.first + 1; rhoIdx++) {
225  for (int jdx = curBin.second - 1; jdx <= curBin.second + 1; jdx++) {
226  // Skip the self reference
227  if (rhoIdx == curBin.first && jdx == curBin.second) continue;
228 
229  // Theta bin needs to handle the wrap.
230  int thetaIdx(jdx);
231 
232  if (thetaIdx < 0)
233  thetaIdx = m_thetaBins - 1;
234  else if (thetaIdx > m_thetaBins - 1)
235  thetaIdx = 0;
236 
237  BinIndex binIndex(rhoIdx, thetaIdx);
238  RhoThetaAccumulatorBinMap::iterator mapItr = rhoThetaAccumulatorBinMap.find(binIndex);
239 
240  if (mapItr != rhoThetaAccumulatorBinMap.end()) {
241  if (mapItr->second.getAccumulatorValues().size() >= threshold)
242  neighborPts.push_back(binIndex);
243  }
244  }
245  }
246 
247  return;
248  }
249 
250  void
252  HoughCluster& neighborPts,
253  HoughCluster& houghCluster,
254  RhoThetaAccumulatorBinMap& rhoThetaAccumulatorBinMap,
255  size_t threshold) const
256  {
257  /**
258  * @brief The workhorse routine for a DBScan like clustering routine to identify peak bins in Hough Space
259  */
260 
261  // Start by adding the input point to our Hough Cluster
262  houghCluster.push_back(curBin);
263 
264  for (auto& binIndex : neighborPts) {
265  AccumulatorBin& accBin = rhoThetaAccumulatorBinMap[binIndex];
266 
267  if (!accBin.isVisited()) {
268  accBin.setVisited();
269 
270  HoughCluster nextNeighborPts;
271 
272  HoughRegionQuery(binIndex, rhoThetaAccumulatorBinMap, nextNeighborPts, threshold);
273 
274  if (!nextNeighborPts.empty()) {
275  for (auto& neighborIdx : nextNeighborPts) {
276  neighborPts.push_back(neighborIdx);
277  }
278  }
279  }
280 
281  if (!accBin.isInCluster()) {
282  houghCluster.push_back(binIndex);
283  accBin.setInCluster();
284  }
285  }
286 
287  return;
288  }
289 
290  bool
292  {
293  return *left < *right;
294  }
295 
296  struct Hit3DSetCompare {
297  bool
298  operator()(const reco::ClusterHit3D* left, const reco::ClusterHit3D* right) const
299  {
300  return Hit3DCompare(left, right);
301  }
302  };
303 
304  class OrderHitsAlongWire {
305  public:
306  OrderHitsAlongWire(int plane = 0) : m_plane(plane) {}
307 
308  bool
310  {
311  int planeToCheck = (m_plane + 1) % 3;
312 
313  return left->getHits()[planeToCheck]->WireID().Wire <
314  right->getHits()[planeToCheck]->WireID().Wire;
315  }
316 
317  private:
318  int m_plane;
319  };
320 
322  bool
323  operator()(const std::pair<size_t, size_t>& left, const std::pair<size_t, size_t>& right)
324  {
325  return left.second < right.second;
326  }
327  };
328 
329  void
331  reco::HitPairListPtr& outputList) const
332  {
333  // Intention is to try to find the largest "contiguous" block of hits in the input list
334  // In a nutshell, the idea is to order input hits according to the pca, then
335  // loop through the hits and store them in a new hit list. If a gap is detected then
336  // we terminate the previous list, create a new one and continue. After the loop over
337  // hits is complete then simply pick the biggest list.
338  // Step I is to run the pca and order the hits
340 
341  m_pcaAlg.PCAAnalysis_3D(inputHitList, pca, true);
342 
343  // It would seem that the analysis can fail!
344  if (!pca.getSvdOK()) {
345  outputList = inputHitList;
346  return;
347  }
348 
349  m_pcaAlg.PCAAnalysis_calc3DDocas(inputHitList, pca);
350 
351  inputHitList.sort(SeedFinderAlgBase::Sort3DHitsByArcLen3D());
352  outputList.clear();
353 
354  // Create containers to hold our hit lists
355  reco::HitPairListPtr continuousHitList;
356 
357  std::map<size_t, reco::HitPairListPtr> gapHitListMap;
358 
359  // Remember the distance arc length of the last hit
360  double arcLenLastHit = inputHitList.front()->getArclenToPoca();
361 
362  // Loop through the input hits
363  for (const auto& hit3D : inputHitList) {
364  // Hits in order, delta arclength should always be positive
365  double arcLen = hit3D->getArclenToPoca();
366  double deltaArcLen = arcLen - arcLenLastHit;
367 
368  if (deltaArcLen > m_maximumGap) {
369  gapHitListMap[continuousHitList.size()] = continuousHitList;
370  continuousHitList.clear();
371  }
372 
373  continuousHitList.emplace_back(hit3D);
374 
375  arcLenLastHit = arcLen;
376  }
377 
378  if (!continuousHitList.empty()) gapHitListMap[continuousHitList.size()] = continuousHitList;
379 
380  // Get the largest list of hits
381  std::map<size_t, reco::HitPairListPtr>::reverse_iterator longestListItr =
382  gapHitListMap.rbegin();
383 
384  if (longestListItr != gapHitListMap.rend()) {
385  size_t nContinuousHits = longestListItr->first;
386  reco::HitPairListPtr longestList = longestListItr->second;
387 
388  outputList.resize(nContinuousHits);
389  std::copy(longestList.begin(), longestList.end(), outputList.begin());
390  }
391 
392  return;
393  }
394 
395  void
398  int& nLoops,
399  RhoThetaAccumulatorBinMap& rhoThetaAccumulatorBinMap,
400  HoughClusterList& houghClusters) const
401  {
402  // The goal of this function is to do a basic Hough Transform on the input list of 3D hits.
403  // In order to transform this to a 2D problem, the 3D hits are projected to the plane of the two
404  // largest eigen values from the input principal components analysis axis.
405  // There are two basic steps to the job here:
406  // 1) Build and accumlate a rho-theta map
407  // 2) Go through that rho-theta map and find candidate Hough "clusters"
408  // Unfortunately, the following may not be suitable viewing for those who may be feint of heart
409  //
410  // Define some constants
411  static int histCount(0);
412  const double maxTheta(M_PI); // Obviously, 180 degrees
413  const double thetaBinSize(maxTheta / double(m_thetaBins)); // around 4 degrees (45 bins)
414  const double rhoBinSizeMin(m_geometry->WirePitch()); // Wire spacing gives a natural bin size?
415 
416  // Recover the parameters from the Principal Components Analysis that we need to project and accumulate
417  Eigen::Vector3f pcaCenter(
418  pca.getAvePosition()[0], pca.getAvePosition()[1], pca.getAvePosition()[2]);
419  Eigen::Vector3f planeVec0(pca.getEigenVectors().row(2));
420  Eigen::Vector3f planeVec1(pca.getEigenVectors().row(1));
421  Eigen::Vector3f pcaPlaneNrml(pca.getEigenVectors().row(0));
422  double eigenVal0 = 3. * sqrt(pca.getEigenValues()[2]);
423  double eigenVal1 = 3. * sqrt(pca.getEigenValues()[1]);
424  double maxRho = std::sqrt(eigenVal0 * eigenVal0 + eigenVal1 * eigenVal1) * 2. / 3.;
425  double rhoBinSize = maxRho / double(m_rhoBins);
426 
427  if (rhoBinSize < rhoBinSizeMin) rhoBinSize = rhoBinSizeMin;
428 
429  // **********************************************************************
430  // Part I: Accumulate values in the rho-theta map
431  // **********************************************************************
432 
433  size_t maxBinCount(0);
434  int nAccepted3DHits(0);
435 
436  // Commence looping over the input 3D hits and fill our accumulator bins
437  for (reco::HitPairListPtr::const_iterator hit3DItr = hitPairListPtr.begin();
438  hit3DItr != hitPairListPtr.end();
439  hit3DItr++) {
440  // Recover the hit
441  const reco::ClusterHit3D* hit3D(*hit3DItr);
442 
443  // Skip hits which are not skeleton points
444  if (!(hit3D->getStatusBits() & 0x10000000)) continue;
445 
446  nAccepted3DHits++;
447 
448  Eigen::Vector3f hit3DPosition(
449  hit3D->getPosition()[0], hit3D->getPosition()[1], hit3D->getPosition()[2]);
450  Eigen::Vector3f pcaToHitVec = hit3DPosition - pcaCenter;
451  Eigen::Vector3f pcaToHitPlaneVec(pcaToHitVec.dot(planeVec0), pcaToHitVec.dot(planeVec1), 0.);
452  double xPcaToHit = pcaToHitPlaneVec[0];
453  double yPcaToHit = pcaToHitPlaneVec[1];
454 
455  // Create an accumulator value
456  AccumulatorValues accValue(pcaToHitPlaneVec, hit3DItr);
457 
458  // Commence loop over theta to fill accumulator bins
459  // Note that with theta in the range 0-pi then we can have negative values for rho
460  for (int thetaIdx = 0; thetaIdx < m_thetaBins; thetaIdx++) {
461  // We need to convert our theta index to an angle
462  double theta = thetaBinSize * double(thetaIdx);
463 
464  // calculate rho for this angle
465  double rho = xPcaToHit * std::cos(theta) + yPcaToHit * std::sin(theta);
466 
467  // Get the rho index
468  int rhoIdx = std::round(rho / rhoBinSize);
469 
470  // Accumulate
471  BinIndex binIndex(rhoIdx, thetaIdx);
472 
473  rhoThetaAccumulatorBinMap[binIndex].addAccumulatorValue(accValue);
474 
475  if (rhoThetaAccumulatorBinMap[binIndex].getAccumulatorValues().size() > maxBinCount)
476  maxBinCount = rhoThetaAccumulatorBinMap[binIndex].getAccumulatorValues().size();
477  }
478  }
479 
480  // Accumulation done, if asked now display the hist
481  if (m_displayHist) {
482  std::ostringstream ostr;
483  ostr << "Hough Histogram " << histCount++;
484  m_Canvases.emplace_back(new TCanvas(ostr.str().c_str(), ostr.str().c_str(), 1000, 1000));
485 
486  std::ostringstream ostr2;
487  ostr2 << "Plane";
488 
489  m_Canvases.back()->GetFrame()->SetFillColor(46);
490  m_Canvases.back()->SetFillColor(19);
491  m_Canvases.back()->SetBorderMode(19);
492  m_Canvases.back()->cd(1);
493 
494  double zmin = 0.06;
495  double zmax = 0.94;
496  double xmin = 0.04;
497  double xmax = 0.95;
498  TPad* p = new TPad(ostr2.str().c_str(), ostr2.str().c_str(), zmin, xmin, zmax, xmax);
499  p->SetBit(kCanDelete); // Give away ownership.
500  p->Range(zmin, xmin, zmax, xmax);
501  p->SetFillStyle(4000); // Transparent.
502  p->Draw();
503  m_Pads.push_back(p);
504 
505  TH2D* houghHist = new TH2D("HoughHist",
506  "Hough Space",
507  2 * m_rhoBins,
508  -m_rhoBins + 0.5,
509  m_rhoBins + 0.5,
510  m_thetaBins,
511  0.,
512  m_thetaBins);
513 
514  for (const auto& rhoThetaMap : rhoThetaAccumulatorBinMap) {
515  houghHist->Fill(rhoThetaMap.first.first,
516  rhoThetaMap.first.second + 0.5,
517  rhoThetaMap.second.getAccumulatorValues().size());
518  }
519 
520  houghHist->SetBit(kCanDelete);
521  houghHist->Draw();
522  m_Canvases.back()->Update();
523  }
524 
525  // **********************************************************************
526  // Part II: Use DBScan (or a slight variation) to find clusters of bins
527  // **********************************************************************
528 
529  size_t thresholdLo = std::max(size_t(m_hiThresholdFrac * nAccepted3DHits), m_hiThresholdMin);
530  size_t thresholdHi = m_loThresholdFrac * maxBinCount;
531 
532  std::list<RhoThetaAccumulatorBinMap::iterator> binIndexList;
533 
534  for (RhoThetaAccumulatorBinMap::iterator mapItr = rhoThetaAccumulatorBinMap.begin();
535  mapItr != rhoThetaAccumulatorBinMap.end();
536  mapItr++)
537  binIndexList.push_back(mapItr);
538 
539  binIndexList.sort(SortBinIndexList());
540 
541  for (auto& mapItr : binIndexList) {
542  // If we have been here before we skip
543  //if (mapItr.second.isVisited()) continue;
544  if (mapItr->second.isInCluster()) continue;
545 
546  // Mark this bin as visited
547  // Actually, don't mark it since we are double thresholding and don't want it missed
548  //mapItr.second.setVisited();
549 
550  // Make sure over threshold
551  if (mapItr->second.getAccumulatorValues().size() < thresholdLo) {
552  mapItr->second.setNoise();
553  continue;
554  }
555 
556  // Set the low threshold to make sure we merge bins that might be either side of a boundary trajectory
557  thresholdHi = std::max(
558  size_t(m_loThresholdFrac * mapItr->second.getAccumulatorValues().size()), m_hiThresholdMin);
559 
560  // Recover our neighborhood
561  HoughCluster neighborhood;
562  BinIndex curBin(mapItr->first);
563 
564  HoughRegionQuery(curBin, rhoThetaAccumulatorBinMap, neighborhood, thresholdHi);
565 
566  houghClusters.push_back(HoughCluster());
567 
568  HoughCluster& houghCluster = houghClusters.back();
569 
571  curBin, neighborhood, houghCluster, rhoThetaAccumulatorBinMap, thresholdHi);
572  }
573 
574  // Sort the clusters using the SortHoughClusterList metric
575  if (!houghClusters.empty()) houghClusters.sort(SortHoughClusterList(rhoThetaAccumulatorBinMap));
576 
577  return;
578  }
579 
580  bool
582  SeedHitPairListPair& seedHitPair) const
583  {
584  if (seed3DHits.size() < m_minimum3DHits) return false;
585 
586  reco::PrincipalComponents seedFullPca;
587 
588  m_pcaAlg.PCAAnalysis_3D(seed3DHits, seedFullPca, true);
589 
590  if (!seedFullPca.getSvdOK()) return false;
591 
592  // Use the following to set the 3D doca and arclength for each hit
593  m_pcaAlg.PCAAnalysis_calc3DDocas(seed3DHits, seedFullPca);
594 
595  // Use this info to sort the hits along the principle axis
596  //seed3DHits.sort(SeedFinderAlgBase::Sort3DHitsByArcLen3D());
598 
599  // The idea here is to search for the first hit that lies "close" to the principle axis
600  // At that point we count out n hits to use as the seed
601  reco::HitPairListPtr seedHit3DList;
602  std::set<const reco::ClusterHit2D*> seedHitSet;
603  double aveDocaToAxis = seedFullPca.getAveHitDoca();
604  int gapCount(0);
605 
606  // Now loop through hits to search for a "continuous" block of at least m_numSeed2DHits
607  // We'll arrive at that number by collecting 2D hits in an stl set which will keep track of unique occurances
608  for (reco::HitPairListPtr::iterator peakBinItr = seed3DHits.begin();
609  peakBinItr != seed3DHits.end();
610  peakBinItr++) {
611  const reco::ClusterHit3D* hit3D = *peakBinItr;
612 
613  if (hit3D->getDocaToAxis() < m_numAveDocas * aveDocaToAxis) {
614  // Check if we need to reset because of gap count
615  if (gapCount > m_numSkippedHits) {
616  seedHit3DList.clear();
617  seedHitSet.clear();
618  }
619 
620  seedHit3DList.push_back(hit3D);
621 
622  for (const auto& hit : hit3D->getHits())
623  seedHitSet.insert(hit);
624 
625  gapCount = 0;
626  }
627  else
628  gapCount++;
629 
630  if (seedHitSet.size() > m_numSeed2DHits) break;
631  }
632 
633  // If not enough hits then we are done
634  if (seedHit3DList.size() < m_minimum3DHits) return false;
635 
637 
638  // Use only the "seed" 3D hits to get new PCA axes
639  m_pcaAlg.PCAAnalysis_3D(seedHit3DList, seedPca, true);
640 
641  if (!seedPca.getSvdOK()) return false;
642 
643  m_pcaAlg.PCAAnalysis_calc3DDocas(seedHit3DList, seedPca);
644  //seedHit3DList.sort(SeedFinderAlgBase::Sort3DHitsByAbsArcLen3D());
645  seedHit3DList.sort(SeedFinderAlgBase::Sort3DHitsByArcLen3D());
646 
647  // Now translate the seedCenter by the arc len to the first hit
648  double seedCenter[3] = {
649  seedPca.getAvePosition()[0], seedPca.getAvePosition()[1], seedPca.getAvePosition()[2]};
650  double seedDir[3] = {seedPca.getEigenVectors().row(2)[0],
651  seedPca.getEigenVectors().row(2)[1],
652  seedPca.getEigenVectors().row(2)[2]};
653 
654  double arcLen = seedHit3DList.front()->getArclenToPoca();
655  double seedStart[3] = {seedCenter[0] + arcLen * seedDir[0],
656  seedCenter[1] + arcLen * seedDir[1],
657  seedCenter[2] + arcLen * seedDir[2]};
658 
659  //seedStart[0] = seedHit3DList.front()->getX();
660  //seedStart[1] = seedHit3DList.front()->getY();
661  //seedStart[2] = seedHit3DList.front()->getZ();
662 
663  if (seedHitSet.size() >= 10) {
664  TVector3 newSeedPos;
665  TVector3 newSeedDir;
666  double chiDOF;
667 
668  LineFit2DHits(seedHitSet, seedStart[0], newSeedPos, newSeedDir, chiDOF);
669 
670  if (chiDOF > 0.) {
671  // check angles between new/old directions
672  double cosAng =
673  seedDir[0] * newSeedDir[0] + seedDir[1] * newSeedDir[1] + seedDir[2] * newSeedDir[2];
674 
675  if (cosAng < 0.) newSeedDir *= -1.;
676 
677  seedStart[0] = newSeedPos[0];
678  seedStart[1] = newSeedPos[1];
679  seedStart[2] = newSeedPos[2];
680  seedDir[0] = newSeedDir[0];
681  seedDir[1] = newSeedDir[1];
682  seedDir[2] = newSeedDir[2];
683  }
684  }
685 
686  // Keep track of this seed and the 3D hits that make it up
687  seedHitPair = SeedHitPairListPair(recob::Seed(seedStart, seedDir), seedHit3DList);
688 
689  // We are going to drop a few hits off the ends in the hope this facilitates finding more track like objects, provided there are enough hits
690  if (seed3DHits.size() > 100) {
691  // Need to reset the doca/arclen first
692  m_pcaAlg.PCAAnalysis_calc3DDocas(seed3DHits, seedFullPca);
693 
694  // Now resort the hits
696 
697  size_t numToKeep = seed3DHits.size() - 10;
698 
699  reco::HitPairListPtr::iterator endItr = seed3DHits.begin();
700 
701  std::advance(endItr, numToKeep);
702 
703  seed3DHits.erase(endItr, seed3DHits.end());
704  }
705 
706  return true;
707  }
708 
709  bool
711  reco::PrincipalComponents& inputPCA,
712  SeedHitPairListPairVec& seedHitPairVec) const
713  {
714  // This will be a busy routine... the basic tasks are:
715  // 1) loop through hits and project to the plane defined by the two largest eigen values, accumulate in Hough space
716  // 2) "Cluster" the Hough space to associate hits which are common to a line
717  // 3) Process these clusters (still to be defined exactly)
718 
719  // Create an interim data structure which will allow us to sort our seeds by "best"
720  // before we return them in the seedHitPairVec
721  typedef std::map<size_t, SeedHitPairListPairVec> SizeToSeedToHitMap;
722 
723  SizeToSeedToHitMap seedHitPairMap;
724  SeedHitPairListPair seedHitPair;
725 
726  // Make sure we are using the right pca
727  reco::HitPairListPtr hitPairListPtr = inputHitPairListPtr;
728 
729  int nLoops(0);
730 
731  // Make a local copy of the input PCA
732  reco::PrincipalComponents pca = inputPCA;
733 
734  // We loop over hits in our list until there are no more
735  while (!hitPairListPtr.empty()) {
736  // We also require that there be some spread in the data, otherwise not worth running?
737  double eigenVal0 = 3. * sqrt(pca.getEigenValues()[2]);
738  double eigenVal1 = 3. * sqrt(pca.getEigenValues()[1]);
739 
740  if (eigenVal0 > 5. && eigenVal1 > 0.001) {
741  // **********************************************************************
742  // Part I: Build Hough space and find Hough clusters
743  // **********************************************************************
744  RhoThetaAccumulatorBinMap rhoThetaAccumulatorBinMap;
745  HoughClusterList houghClusters;
746 
747  findHoughClusters(hitPairListPtr, pca, nLoops, rhoThetaAccumulatorBinMap, houghClusters);
748 
749  // If no clusters then done
750  if (houghClusters.empty()) break;
751 
752  // **********************************************************************
753  // Part II: Go through the clusters to find the peak bins
754  // **********************************************************************
755 
756  // We need to use a set so we can be sure to have unique hits
757  reco::HitPairListPtr clusterHitsList;
758  std::set<const reco::ClusterHit3D*> masterHitPtrList;
759  std::set<const reco::ClusterHit3D*> peakBinPtrList;
760 
761  size_t firstPeakCount(0);
762 
763  // Loop through the list of all clusters found above
764  for (auto& houghCluster : houghClusters) {
765  BinIndex peakBin = houghCluster.front();
766  size_t peakCount = 0;
767  size_t totalHits = 0;
768 
769  // Make a local (to this cluster) set of of hits
770  std::set<const reco::ClusterHit3D*> localHitPtrList;
771 
772  // Now loop through the bins that were attached to this cluster
773  for (auto& binIndex : houghCluster) {
774  // An even more local list so we can keep track of peak values
775  std::set<const reco::ClusterHit3D*> tempHitPtrList;
776 
777  // Recover the hits associated to this cluster
778  for (auto& hitItr : rhoThetaAccumulatorBinMap[binIndex].getAccumulatorValues()) {
779  reco::HitPairListPtr::const_iterator hit3DItr = hitItr.getHitIterator();
780 
781  tempHitPtrList.insert(*hit3DItr);
782  }
783 
784  // count hits before we remove any
785  totalHits += tempHitPtrList.size();
786 
787  // Trim out any hits already used by a bigger/better cluster
788  std::set<const reco::ClusterHit3D*> tempHit3DSet;
789 
790  std::set_difference(tempHitPtrList.begin(),
791  tempHitPtrList.end(),
792  masterHitPtrList.begin(),
793  masterHitPtrList.end(),
794  std::inserter(tempHit3DSet, tempHit3DSet.end()));
795 
796  tempHitPtrList = tempHit3DSet;
797 
798  size_t binCount = tempHitPtrList.size();
799 
800  if (peakCount < binCount) {
801  peakCount = binCount;
802  peakBin = binIndex;
803  peakBinPtrList = tempHitPtrList;
804  }
805 
806  // Add this to our local list
807  localHitPtrList.insert(tempHitPtrList.begin(), tempHitPtrList.end());
808  }
809 
810  if (localHitPtrList.size() < m_minimum3DHits) continue;
811 
812  if (!firstPeakCount) firstPeakCount = peakCount;
813 
814  // If the peak counts are significantly less than the first cluster's peak then skip
815  if (peakCount < firstPeakCount / 10) continue;
816 
817  // **********************************************************************
818  // Part III: Make a Seed from the peak bin hits
819  // **********************************************************************
820 
821  reco::HitPairListPtr allPeakBinHits;
822 
823  for (const auto& hit3D : localHitPtrList)
824  allPeakBinHits.push_back(hit3D);
825 
826  reco::HitPairListPtr peakBinHits;
827 
828  // Find longest "continuous" set of hits and use these for the seed
829  findHitGaps(allPeakBinHits, peakBinHits);
830 
831  if (peakBinHits.size() < m_minimum3DHits) continue;
832 
833  // We now build the actual seed.
834  if (buildSeed(peakBinHits, seedHitPair)) {
835  // Keep track of this in our map (which will do ordering for us)
836  seedHitPairMap[peakBinHits.size()].push_back(seedHitPair);
837 
838  // For visual testing in event display, mark all the hits in the first seed so we can see them
839  if (seedHitPairMap.size() == 1) {
840  for (const auto& hit3D : peakBinHits)
841  hit3D->setStatusBit(0x40000000);
842  }
843 
844  // for(const auto& hit3D : seedHitPair.second) hit3D->setStatusBit(0x40000000);
845  }
846 
847  // Our peakBinHits collection will most likely be a subset of the localHitPtrList collection
848  // We want to remove only the "pure" hits which are those in the peakBinHits collection
849  // So, sort them and then add to our master list
850  peakBinHits.sort();
851 
852  masterHitPtrList.insert(peakBinHits.begin(), peakBinHits.end());
853 
854  if (hitPairListPtr.size() - masterHitPtrList.size() < m_minimum3DHits) break;
855  } // end loop over hough clusters
856 
857  // If the masterHitPtrList is empty then nothing happened and we're done
858  if (masterHitPtrList.empty()) break;
859 
860  // **********************************************************************
861  // Part IV: Remove remaining peak bin hits from HitPairPtrList
862  // **********************************************************************
863 
864  hitPairListPtr.sort();
865 
866  reco::HitPairListPtr::iterator newListEnd = std::set_difference(hitPairListPtr.begin(),
867  hitPairListPtr.end(),
868  masterHitPtrList.begin(),
869  masterHitPtrList.end(),
870  hitPairListPtr.begin());
871 
872  hitPairListPtr.erase(newListEnd, hitPairListPtr.end());
873 
874  if (hitPairListPtr.size() < m_minimum3DHits) break;
875 
876  if (nLoops++ > m_maxLoopsPerCluster) break;
877 
878  // ********************************************************
879  }
880  else
881  break; // eigen values not in range
882 
883  // At this point run the PCA on the remaining hits
884  m_pcaAlg.PCAAnalysis_3D(hitPairListPtr, pca, true);
885 
886  if (!pca.getSvdOK()) break;
887  }
888 
889  // The final task before returning is to transfer the stored seeds into the output seed vector
890  // What we want to do is make sure the first seeds are the "best" seeds which is defined as the
891  // seeds which were associated to the most hits by the Hough Transform. Our seed map will have
892  // the reverse of this ordering so we simply iterate through it "backwards"
893  for (SizeToSeedToHitMap::reverse_iterator seedMapItr = seedHitPairMap.rbegin();
894  seedMapItr != seedHitPairMap.rend();
895  seedMapItr++) {
896  for (const auto& seedHitPair : seedMapItr->second) {
897  seedHitPairVec.emplace_back(seedHitPair);
898  }
899  }
900 
901  return true;
902  }
903 
904  bool
906  reco::PrincipalComponents& inputPCA,
907  reco::HitPairListPtrList& hitPairListPtrList) const
908  {
909  // The goal of this routine is run the Hough Transform on the input set of hits
910  // and then to return a list of lists of hits which are associated to a given line
911 
912  // Make sure we are using the right pca
913  reco::HitPairListPtr hitPairListPtr = inputHitPairListPtr;
914 
915  int nLoops(0);
916 
917  // Make a local copy of the input PCA
918  reco::PrincipalComponents pca = inputPCA;
919 
920  // We also require that there be some spread in the data, otherwise not worth running?
921  double eigenVal0 = 3. * sqrt(pca.getEigenValues()[2]);
922  double eigenVal1 = 3. * sqrt(pca.getEigenValues()[1]);
923 
924  if (eigenVal0 > 5. && eigenVal1 > 0.001) {
925  // **********************************************************************
926  // Part I: Build Hough space and find Hough clusters
927  // **********************************************************************
928  RhoThetaAccumulatorBinMap rhoThetaAccumulatorBinMap;
929  HoughClusterList houghClusters;
930 
931  findHoughClusters(hitPairListPtr, pca, nLoops, rhoThetaAccumulatorBinMap, houghClusters);
932 
933  // **********************************************************************
934  // Part II: Go through the clusters to find the peak bins
935  // **********************************************************************
936 
937  // We need to use a set so we can be sure to have unique hits
938  reco::HitPairListPtr clusterHitsList;
939  std::set<const reco::ClusterHit3D*> masterHitPtrList;
940  std::set<const reco::ClusterHit3D*> peakBinPtrList;
941 
942  size_t firstPeakCount(0);
943 
944  // Loop through the list of all clusters found above
945  for (auto& houghCluster : houghClusters) {
946  BinIndex peakBin = houghCluster.front();
947  size_t peakCount = 0;
948  size_t totalHits = 0;
949 
950  // Make a local (to this cluster) set of of hits
951  std::set<const reco::ClusterHit3D*> localHitPtrList;
952 
953  // Now loop through the bins that were attached to this cluster
954  for (auto& binIndex : houghCluster) {
955  // An even more local list so we can keep track of peak values
956  std::set<const reco::ClusterHit3D*> tempHitPtrList;
957 
958  // Recover the hits associated to this cluster
959  for (auto& hitItr : rhoThetaAccumulatorBinMap[binIndex].getAccumulatorValues()) {
960  reco::HitPairListPtr::const_iterator hit3DItr = hitItr.getHitIterator();
961 
962  tempHitPtrList.insert(*hit3DItr);
963  }
964 
965  // count hits before we remove any
966  totalHits += tempHitPtrList.size();
967 
968  // Trim out any hits already used by a bigger/better cluster
969  std::set<const reco::ClusterHit3D*> tempHit3DSet;
970 
971  std::set_difference(tempHitPtrList.begin(),
972  tempHitPtrList.end(),
973  masterHitPtrList.begin(),
974  masterHitPtrList.end(),
975  std::inserter(tempHit3DSet, tempHit3DSet.end()));
976 
977  tempHitPtrList = tempHit3DSet;
978 
979  size_t binCount = tempHitPtrList.size();
980 
981  if (peakCount < binCount) {
982  peakCount = binCount;
983  peakBin = binIndex;
984  peakBinPtrList = tempHitPtrList;
985  }
986 
987  // Add this to our local list
988  localHitPtrList.insert(tempHitPtrList.begin(), tempHitPtrList.end());
989  }
990 
991  if (localHitPtrList.size() < m_minimum3DHits) continue;
992 
993  if (!firstPeakCount) firstPeakCount = peakCount;
994 
995  // If the peak counts are significantly less than the first cluster's peak then skip
996  if (peakCount < firstPeakCount / 10) continue;
997 
998  // **********************************************************************
999  // Part III: Make a list of hits from the total number associated
1000  // **********************************************************************
1001 
1002  hitPairListPtrList.push_back(reco::HitPairListPtr());
1003 
1004  hitPairListPtrList.back().resize(localHitPtrList.size());
1005  std::copy(
1006  localHitPtrList.begin(), localHitPtrList.end(), hitPairListPtrList.back().begin());
1007 
1008  // We want to remove the hits which have been used from further contention
1009  masterHitPtrList.insert(localHitPtrList.begin(), localHitPtrList.end());
1010 
1011  if (hitPairListPtr.size() - masterHitPtrList.size() < m_minimum3DHits) break;
1012  } // end loop over hough clusters
1013  }
1014 
1015  return true;
1016  }
1017 
1018  //------------------------------------------------------------------------------
1019  void
1020  HoughSeedFinderAlg::LineFit2DHits(std::set<const reco::ClusterHit2D*>& hit2DSet,
1021  double XOrigin,
1022  TVector3& Pos,
1023  TVector3& Dir,
1024  double& ChiDOF) const
1025  {
1026  // The following is lifted from Bruce Baller to try to get better
1027  // initial parameters for a candidate Seed. It is slightly reworked
1028  // which is why it is included here instead of used as is.
1029  //
1030  // Linear fit using X as the independent variable. Hits to be fitted
1031  // are passed in the hits vector in a pair form (X, WireID). The
1032  // fitted track position at XOrigin is returned in the Pos vector.
1033  // The direction cosines are returned in the Dir vector.
1034  //
1035  // SVD fit adapted from $ROOTSYS/tutorials/matrix/solveLinear.C
1036  // Fit equation is w = A(X)v, where w is a vector of hit wires, A is
1037  // a matrix to calculate a track projected to a point at X, and v is
1038  // a vector (Yo, Zo, dY/dX, dZ/dX).
1039  //
1040  // Note: The covariance matrix should also be returned
1041  // B. Baller August 2014
1042 
1043  // assume failure
1044  ChiDOF = -1;
1045 
1046  if (hit2DSet.size() < 4) return;
1047 
1048  const unsigned int nvars = 4;
1049  unsigned int npts = hit2DSet.size();
1050 
1051  TMatrixD A(npts, nvars);
1052  // vector holding the Wire number
1053  TVectorD w(npts);
1054  unsigned short ninpl[3] = {0};
1055  unsigned short nok = 0;
1056  unsigned short iht(0), cstat, tpc, ipl;
1057  double x, cw, sw, off;
1058 
1059  // Loop over unique 2D hits from the input list of 3D hits
1060  for (const auto& hit : hit2DSet) {
1061  geo::WireID wireID = hit->WireID();
1062 
1063  cstat = wireID.Cryostat;
1064  tpc = wireID.TPC;
1065  ipl = wireID.Plane;
1066 
1067  // get the wire plane offset
1068  off = m_geometry->WireCoordinate(0, 0, ipl, tpc, cstat);
1069 
1070  // get the "cosine-like" component
1071  cw = m_geometry->WireCoordinate(1, 0, ipl, tpc, cstat) - off;
1072 
1073  // the "sine-like" component
1074  sw = m_geometry->WireCoordinate(0, 1, ipl, tpc, cstat) - off;
1075 
1076  x = hit->getXPosition() - XOrigin;
1077 
1078  A[iht][0] = cw;
1079  A[iht][1] = sw;
1080  A[iht][2] = cw * x;
1081  A[iht][3] = sw * x;
1082  w[iht] = wireID.Wire - off;
1083 
1084  ++ninpl[ipl];
1085 
1086  // need at least two points in a plane
1087  if (ninpl[ipl] == 2) ++nok;
1088 
1089  iht++;
1090  }
1091 
1092  // need at least 2 planes with at least two points
1093  if (nok < 2) return;
1094 
1095  TDecompSVD svd(A);
1096  bool ok;
1097  TVectorD tVec = svd.Solve(w, ok);
1098 
1099  ChiDOF = 0;
1100 
1101  // not enough points to calculate Chisq
1102  if (npts <= 4) return;
1103 
1104  double ypr, zpr, diff;
1105 
1106  for (const auto& hit : hit2DSet) {
1107  geo::WireID wireID = hit->WireID();
1108 
1109  cstat = wireID.Cryostat;
1110  tpc = wireID.TPC;
1111  ipl = wireID.Plane;
1112  off = m_geometry->WireCoordinate(0, 0, ipl, tpc, cstat);
1113  cw = m_geometry->WireCoordinate(1, 0, ipl, tpc, cstat) - off;
1114  sw = m_geometry->WireCoordinate(0, 1, ipl, tpc, cstat) - off;
1115  x = hit->getXPosition() - XOrigin;
1116  ypr = tVec[0] + tVec[2] * x;
1117  zpr = tVec[1] + tVec[3] * x;
1118  diff = ypr * cw + zpr * sw - (wireID.Wire - off);
1119  ChiDOF += diff * diff;
1120  }
1121 
1122  float werr2 = m_geometry->WirePitch() * m_geometry->WirePitch();
1123  ChiDOF /= werr2;
1124  ChiDOF /= (float)(npts - 4);
1125 
1126  double norm = sqrt(1 + tVec[2] * tVec[2] + tVec[3] * tVec[3]);
1127  Dir[0] = 1 / norm;
1128  Dir[1] = tVec[2] / norm;
1129  Dir[2] = tVec[3] / norm;
1130 
1131  Pos[0] = XOrigin;
1132  Pos[1] = tVec[0];
1133  Pos[2] = tVec[1];
1134 
1135  } // TrkLineFit()
1136 
1137 } // namespace lar_cluster3d
geo::Length_t WireCoordinate(double YPos, double ZPos, geo::PlaneID const &planeid) const
Returns the index of the nearest wire to the specified position.
void LineFit2DHits(std::set< const reco::ClusterHit2D * > &hitList, double XOrigin, TVector3 &Pos, TVector3 &Dir, double &ChiDOF) const
intermediate_table::iterator iterator
Define a comparator which will sort hits by arc length along a PCA axis.
void PCAAnalysis_calc3DDocas(const reco::HitPairListPtr &hitPairVector, const reco::PrincipalComponents &pca) const
void findHitGaps(reco::HitPairListPtr &inputHitList, reco::HitPairListPtr &outputList) const
Using Principal Components Axis, look for gaps in a list of 3D hits.
bool getSvdOK() const
Definition: Cluster3D.h:244
void PCAAnalysis_3D(const reco::HitPairListPtr &hitPairList, reco::PrincipalComponents &pca, bool skeletonOnly=false) const
T * get() const
Definition: ServiceHandle.h:63
const Eigen::Vector3f getPosition() const
Definition: Cluster3D.h:158
CryostatID_t Cryostat
Index of cryostat.
Definition: geo_types.h:212
std::vector< AccumulatorValues > AccumulatorValuesVec
HoughSeedFinderAlg(fhicl::ParameterSet const &pset)
Constructor.
WireID_t Wire
Index of the wire within its plane.
Definition: geo_types.h:580
intermediate_table::const_iterator const_iterator
std::list< HitPairListPtr > HitPairListPtrList
Definition: Cluster3D.h:337
std::vector< TVirtualPad * > m_Pads
View pads in current canvas.
bool operator()(const reco::ClusterHit3D *left, const reco::ClusterHit3D *right) const
void findHoughClusters(const reco::HitPairListPtr &inputHits, reco::PrincipalComponents &pca, int &nLoops, RhoThetaAccumulatorBinMap &rhoThetaMap, HoughClusterList &clusterList) const
unsigned int getStatusBits() const
Definition: Cluster3D.h:157
geo::Length_t WirePitch(geo::PlaneID const &planeid) const
Returns the distance between two consecutive wires.
std::pair< recob::Seed, reco::HitPairListPtr > SeedHitPairListPair
float getDocaToAxis() const
Definition: Cluster3D.h:169
art framework interface to geometry description
Define a comparator which will sort hits by the absolute value of arc length so hits are ordered clos...
decltype(auto) constexpr size(T &&obj)
ADL-aware version of std::size.
Definition: StdUtils.h:92
std::vector< std::unique_ptr< TCanvas > > m_Canvases
Graphical trace canvases.
void expandHoughCluster(BinIndex &curBin, HoughCluster &neighborPts, HoughCluster &houghCluster, RhoThetaAccumulatorBinMap &rhoThetaAccumulatorBinMap, size_t threshold) const
const EigenValues & getEigenValues() const
Definition: Cluster3D.h:246
void HoughRegionQuery(BinIndex &curBin, RhoThetaAccumulatorBinMap &rhoThetaAccumulatorBinMap, HoughCluster &neighborPts, size_t threshold) const
T get(std::string const &key) const
Definition: ParameterSet.h:271
const Eigen::Vector3f & getAvePosition() const
Definition: Cluster3D.h:248
std::list< HoughCluster > HoughClusterList
p
Definition: test.py:223
const AccumulatorValuesVec & getAccumulatorValues() const
std::list< const reco::ClusterHit3D * > HitPairListPtr
Definition: Cluster3D.h:335
bool operator()(const HoughSeedFinderAlg::HoughCluster &left, const HoughSeedFinderAlg::HoughCluster &right)
bool operator()(const reco::ClusterHit3D *left, const reco::ClusterHit3D *right)
Q_UINT16 values[128]
static int max(int a, int b)
PlaneID_t Plane
Index of the plane within its TPC.
Definition: geo_types.h:493
bool buildSeed(reco::HitPairListPtr &seed3DHits, SeedHitPairListPair &seedHitPair) const
Given a list of candidate "seed" 3D hits, build the seed and get associated unique 2D hits...
bool operator()(const HoughSeedFinderAlg::RhoThetaAccumulatorBinMap::iterator &left, const HoughSeedFinderAlg::RhoThetaAccumulatorBinMap::iterator &right)
This is used to sort bins in Hough Space.
std::map< BinIndex, AccumulatorBin > RhoThetaAccumulatorBinMap
auto norm(Vector const &v)
Return norm of the specified vector.
Detector simulation of raw signals on wires.
#define M_PI
Definition: includeROOT.h:54
virtual bool findTrackHits(reco::HitPairListPtr &hitPairListPtr, reco::PrincipalComponents &inputPCA, reco::HitPairListPtrList &hitPairListPtrList) const
Given the list of hits this will return the sets of hits which belong on the same line...
const Eigen::Vector3f & getPosition() const
bool operator()(const std::pair< size_t, size_t > &left, const std::pair< size_t, size_t > &right)
const float getAveHitDoca() const
Definition: Cluster3D.h:249
Eigen::Vector3f m_position
We really only need the x,y coordinates here but keep all three for now.
T copy(T const &v)
#define A
Definition: memgrp.cpp:38
list x
Definition: train.py:276
virtual bool findTrackSeeds(reco::HitPairListPtr &hitPairListPtr, reco::PrincipalComponents &inputPCA, SeedHitPairListPairVec &seedHitPairVec) const
Given the list of hits this will search for candidate Seed objects and return them.
SortHoughClusterList(HoughSeedFinderAlg::RhoThetaAccumulatorBinMap &accMap)
This is used to sort "Hough Clusters" by the maximum entries in a bin.
std::vector< SeedHitPairListPair > SeedHitPairListPairVec
auto const & get(AssnsNode< L, R, D > const &r)
Definition: AssnsNode.h:115
TPCID_t TPC
Index of the TPC within its cryostat.
Definition: geo_types.h:406
AccumulatorValues()
A utility class to contain the values of a given "bin" in Hough Space.
AccumulatorBin()
A utility class used to accumulate the above values.
const ClusterHit2DVec & getHits() const
Definition: Cluster3D.h:171
AccumulatorValues(const Eigen::Vector3f &position, const reco::HitPairListPtr::const_iterator &itr)
bool Hit3DCompare(const reco::ClusterHit3D *left, const reco::ClusterHit3D *right)
reco::HitPairListPtr::const_iterator getHitIterator() const
const EigenVectors & getEigenVectors() const
Definition: Cluster3D.h:247