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