LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
lar_cluster3d::ParallelHitsSeedFinderAlg Class Referencefinal

ParallelHitsSeedFinderAlg class. More...

#include "ParallelHitsSeedFinderAlg.h"

Inheritance diagram for lar_cluster3d::ParallelHitsSeedFinderAlg:
lar_cluster3d::SeedFinderAlgBase

Public Member Functions

 ParallelHitsSeedFinderAlg (fhicl::ParameterSet const &pset)
 Constructor. More...
 
bool findTrackSeeds (reco::HitPairListPtr &hitPairListPtr, reco::PrincipalComponents &inputPCA, SeedHitPairListPairVec &seedHitMap) const override
 Given the list of hits this will search for candidate Seed objects and return them. More...
 

Private Attributes

size_t m_maxNumEdgeHits
 Maximum number hits each end of PCA axis. More...
 
double m_gapDistance
 Maximum allowed distance between hits. More...
 
size_t m_numSeed2DHits
 Number 2D seed hits desired. More...
 
PrincipalComponentsAlg m_pcaAlg
 

Detailed Description

ParallelHitsSeedFinderAlg class.

Definition at line 27 of file ParallelHitsSeedFinderAlg.h.

Constructor & Destructor Documentation

lar_cluster3d::ParallelHitsSeedFinderAlg::ParallelHitsSeedFinderAlg ( fhicl::ParameterSet const &  pset)
explicit

Constructor.

Parameters
pset

Definition at line 27 of file ParallelHitsSeedFinderAlg.cxx.

References fhicl::ParameterSet::get(), m_gapDistance, m_maxNumEdgeHits, and m_numSeed2DHits.

28  : m_maxNumEdgeHits(1000)
29  , m_gapDistance(20.)
30  , m_numSeed2DHits(80)
31  , m_pcaAlg(pset.get<fhicl::ParameterSet>("PrincipalComponentsAlg"))
32  {
33  m_maxNumEdgeHits = pset.get<size_t>("MaxNumEdgeHits", 1000);
34  m_gapDistance = pset.get<double>("GapDistance", 20.);
35  m_numSeed2DHits = pset.get<size_t>("NumSeed2DHits", 80);
36  }
double m_gapDistance
Maximum allowed distance between hits.
size_t m_maxNumEdgeHits
Maximum number hits each end of PCA axis.
size_t m_numSeed2DHits
Number 2D seed hits desired.

Member Function Documentation

bool lar_cluster3d::ParallelHitsSeedFinderAlg::findTrackSeeds ( reco::HitPairListPtr hitPairListPtr,
reco::PrincipalComponents inputPCA,
SeedHitPairListPairVec seedHitMap 
) const
overridevirtual

Given the list of hits this will search for candidate Seed objects and return them.

Implements lar_cluster3d::SeedFinderAlgBase.

Definition at line 38 of file ParallelHitsSeedFinderAlg.cxx.

References reco::PrincipalComponents::flipAxis(), reco::PrincipalComponents::getEigenVectors(), reco::PrincipalComponents::getSvdOK(), m_gapDistance, m_maxNumEdgeHits, m_numSeed2DHits, m_pcaAlg, lar_cluster3d::PrincipalComponentsAlg::PCAAnalysis_3D(), and lar_cluster3d::PrincipalComponentsAlg::PCAAnalysis_calc3DDocas().

Referenced by lar_cluster3d::Cluster3D::findTrackSeeds().

41  {
42  // This routine can fail...
43  bool foundGoodSeeds(false);
44 
45  // Make sure we are using the right pca
46  reco::PrincipalComponents pca = inputPCA;
47 
48  if (pca.getSvdOK()) {
49  // Presume CR muons will be "downward going"...
50  if (pca.getEigenVectors().row(2)(1) > 0.) pca.flipAxis(0);
51 
52  // This routine is typically called when there are LOTS of hits... so we are going to try
53  // to reduce the number of operations on the full list as much as possible. However, we
54  // really need the input hit list to be sorted by th input PCA so do that now
55  m_pcaAlg.PCAAnalysis_calc3DDocas(inputHitPairListPtr, pca);
56 
57  // Use this info to sort the hits along the principle axis
58  // Note that this will sort hits from the "middle" to the "outside"
59  inputHitPairListPtr.sort(SeedFinderAlgBase::Sort3DHitsByArcLen3D());
60 
61  // We need to determine the number of hits to drop... and this is really driven by the number of unique
62  // 2D hits we want to keep at one end of the cluster. So, make some containers and do some looping
63  reco::HitPairListPtr seedHit3DList;
64  std::set<const reco::ClusterHit2D*> hit2DSet;
65 
66  // Look for numSeed2DHits which are continuous
67  double lastArcLen = inputHitPairListPtr.front()->getArclenToPoca();
68 
69  for (const auto& hit3D : inputHitPairListPtr) {
70  double arcLen = hit3D->getArclenToPoca();
71 
72  if (fabs(arcLen - lastArcLen) > m_gapDistance) {
73  seedHit3DList.clear();
74  hit2DSet.clear();
75  }
76 
77  seedHit3DList.push_back(hit3D);
78 
79  for (const auto& hit2D : hit3D->getHits()) {
80  hit2DSet.insert(hit2D);
81  }
82 
83  if (hit2DSet.size() > m_numSeed2DHits) break;
84 
85  lastArcLen = arcLen;
86  }
87 
88  // We require a minimum number of seed hits in order to proceed
89  if (hit2DSet.size() > m_numSeed2DHits) {
90  // The above derivation determines the number of hits to keep each side of the cloud
91  size_t num3DHitsToKeep = std::min(2 * seedHit3DList.size(), inputHitPairListPtr.size());
92  size_t numEdgeHits = std::min(size_t(num3DHitsToKeep / 2), m_maxNumEdgeHits);
93 
94  // Get an iterator for the hits
95  reco::HitPairListPtr::iterator edgeHitItr = inputHitPairListPtr.begin();
96 
97  std::advance(edgeHitItr, numEdgeHits);
98 
99  // Make a container for copying in the edge hits and size it to hold all the hits
100  reco::HitPairListPtr hit3DList;
101 
102  hit3DList.resize(2 * numEdgeHits);
103 
104  // Copy the low edge hit range into the list
105  reco::HitPairListPtr::iterator nextHit3DItr =
106  std::copy(inputHitPairListPtr.begin(), edgeHitItr, hit3DList.begin());
107 
108  // reset the "seed hits"
109  seedHit3DList.clear();
110  seedHit3DList.resize(numEdgeHits);
111 
112  std::copy(inputHitPairListPtr.begin(), edgeHitItr, seedHit3DList.begin());
113 
114  // Now advance the iterator into the main container and copy the rest of the elements
115  std::advance(edgeHitItr, inputHitPairListPtr.size() - 2 * numEdgeHits);
116 
117  std::copy(edgeHitItr, inputHitPairListPtr.end(), nextHit3DItr);
118 
119  // At this point we should now have trimmed out the bulk of the 3D hits from the middle
120  // of the input cloud of hits. Next step is to rerun the PCA on our reduced set of hits
122 
123  m_pcaAlg.PCAAnalysis_3D(hit3DList, seedPCA, false);
124 
125  if (seedPCA.getSvdOK()) {
126  // Still looking to point "down"
127  if (seedPCA.getEigenVectors().row(2)(1) > 0.) seedPCA.flipAxis(0);
128 
129  // Recompute the 3D docas and arc lengths
130  m_pcaAlg.PCAAnalysis_calc3DDocas(hit3DList, seedPCA);
131 
132  // Now sort hits in regular order
133  //hit3DList.sort(SeedFinderAlgBase::Sort3DHitsByArcLen3D());
134  seedHit3DList.sort(SeedFinderAlgBase::Sort3DHitsByArcLen3D());
135 
136  // Now translate the seedCenter by the arc len to the first hit
137  double seedDir[3] = {seedPCA.getEigenVectors().row(2)(0),
138  seedPCA.getEigenVectors().row(2)(1),
139  seedPCA.getEigenVectors().row(2)(2)};
140  double seedStart[3] = {seedHit3DList.front()->getX(),
141  seedHit3DList.front()->getY(),
142  seedHit3DList.front()->getZ()};
143 
144  // Move from the first hit to the halfway point but from the first hit, not from the axis position
145  double halfArcLen = 0.5 * fabs(seedHit3DList.back()->getArclenToPoca() -
146  seedHit3DList.front()->getArclenToPoca());
147 
148  seedStart[0] += halfArcLen * seedDir[0];
149  seedStart[1] += halfArcLen * seedDir[1];
150  seedStart[2] += halfArcLen * seedDir[2];
151 
152  for (const auto& hit3D : seedHit3DList)
153  hit3D->setStatusBit(0x40000000);
154 
155  seedHitPairVec.emplace_back(std::pair<recob::Seed, reco::HitPairListPtr>(
156  recob::Seed(seedStart, seedDir), seedHit3DList));
157 
158  inputPCA = seedPCA;
159 
160  foundGoodSeeds = true;
161  }
162  }
163  }
164 
165  return foundGoodSeeds;
166  }
intermediate_table::iterator iterator
void PCAAnalysis_calc3DDocas(const reco::HitPairListPtr &hitPairVector, const reco::PrincipalComponents &pca) const
bool getSvdOK() const
Definition: Cluster3D.h:240
void PCAAnalysis_3D(const reco::HitPairListPtr &hitPairList, reco::PrincipalComponents &pca, bool skeletonOnly=false) const
double m_gapDistance
Maximum allowed distance between hits.
void flipAxis(size_t axis)
Definition: Cluster3D.cxx:232
size_t m_maxNumEdgeHits
Maximum number hits each end of PCA axis.
std::list< const reco::ClusterHit3D * > HitPairListPtr
Definition: Cluster3D.h:326
size_t m_numSeed2DHits
Number 2D seed hits desired.
const EigenVectors & getEigenVectors() const
Definition: Cluster3D.h:243

Member Data Documentation

double lar_cluster3d::ParallelHitsSeedFinderAlg::m_gapDistance
private

Maximum allowed distance between hits.

Definition at line 45 of file ParallelHitsSeedFinderAlg.h.

Referenced by findTrackSeeds(), and ParallelHitsSeedFinderAlg().

size_t lar_cluster3d::ParallelHitsSeedFinderAlg::m_maxNumEdgeHits
private

Maximum number hits each end of PCA axis.

Definition at line 44 of file ParallelHitsSeedFinderAlg.h.

Referenced by findTrackSeeds(), and ParallelHitsSeedFinderAlg().

size_t lar_cluster3d::ParallelHitsSeedFinderAlg::m_numSeed2DHits
private

Number 2D seed hits desired.

Definition at line 46 of file ParallelHitsSeedFinderAlg.h.

Referenced by findTrackSeeds(), and ParallelHitsSeedFinderAlg().

PrincipalComponentsAlg lar_cluster3d::ParallelHitsSeedFinderAlg::m_pcaAlg
private

Definition at line 48 of file ParallelHitsSeedFinderAlg.h.

Referenced by findTrackSeeds().


The documentation for this class was generated from the following files: