PCASeedFinderAlg.cxx
Go to the documentation of this file.
1 /**
2  * @file PCASeedFinderAlg.cxx
3  *
4  * @brief Implementation of the Seed Finder Algorithm
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
8  */
9 
11 
12 // Framework Includes
14 #include "fhiclcpp/ParameterSet.h"
15 
16 // LArSoft includes
19 
20 // ROOT includes
21 #include "TDecompSVD.h"
22 #include "TMatrixD.h"
23 #include "TVector3.h"
24 #include "TVectorD.h"
25 
26 // std includes
27 #include <memory>
28 
29 //------------------------------------------------------------------------------------------------------------------------------------------
30 
31 //------------------------------------------------------------------------------------------------------------------------------------------
32 // implementation follows
33 
34 namespace lar_cluster3d {
35 
37  : m_gapDistance(5.)
38  , m_numSeed2DHits(80)
39  , m_minAllowedCosAng(0.7)
40  , m_pcaAlg(pset.get<fhicl::ParameterSet>("PrincipalComponentsAlg"))
41  {
42  m_gapDistance = pset.get<double>("GapDistance", 5.);
43  m_numSeed2DHits = pset.get<size_t>("NumSeed2DHits", 80);
44  m_minAllowedCosAng = pset.get<double>("MinAllowedCosAng", 0.7);
46  }
47 
48  bool
50  reco::PrincipalComponents& inputPCA,
51  SeedHitPairListPairVec& seedHitPairVec) const
52  {
53  bool foundGoodSeed(false);
54 
55  // Make sure we are using the right pca
56  reco::HitPairListPtr hitPairListPtr = inputHitPairListPtr;
57 
58  // Make a local copy of the input pca
59  reco::PrincipalComponents pca = inputPCA;
60 
61  // We also require that there be some spread in the data, otherwise not worth running?
62  double eigenVal0 = 3. * sqrt(pca.getEigenValues()[2]);
63  double eigenVal1 = 3. * sqrt(pca.getEigenValues()[1]);
64 
65  if (eigenVal0 > 5. && eigenVal1 > 0.001) {
66  // Presume CR muons will be "downward going"...
67  if (pca.getEigenVectors().row(2)(1) > 0.) pca.flipAxis(0);
68 
69  // Use the following to set the 3D doca and arclength for each hit
70  m_pcaAlg.PCAAnalysis_calc3DDocas(hitPairListPtr, pca);
71 
72  // Use this info to sort the hits along the principle axis
73  hitPairListPtr.sort(SeedFinderAlgBase::Sort3DHitsByArcLen3D());
74  //hitPairListPtr.sort(PcaSort3DHitsByAbsArcLen3D());
75 
76  // Make a local copy of the 3D hits
77  reco::HitPairListPtr hit3DList;
78 
79  hit3DList.resize(hitPairListPtr.size());
80 
81  std::copy(hitPairListPtr.begin(), hitPairListPtr.end(), hit3DList.begin());
82 
83  reco::PrincipalComponents seedPca = pca;
84 
85  if (getHitsAtEnd(hit3DList, seedPca)) {
86  // We can use the same routine to check hits at the opposite end to make sure
87  // we have consistency between both ends of the track.
88  // So... follow the same general program as above
89  reco::HitPairListPtr checkList;
90 
91  checkList.resize(hitPairListPtr.size());
92 
93  std::copy(hitPairListPtr.begin(), hitPairListPtr.end(), checkList.begin());
94 
95  std::reverse(checkList.begin(), checkList.end());
96 
97  reco::PrincipalComponents checkPca = pca;
98 
99  if (getHitsAtEnd(checkList, checkPca)) {
100  // We want our seed to be in the middle of our seed points
101  // First step to getting there is to reorder the hits along the
102  // new pca axis
103  m_pcaAlg.PCAAnalysis_calc3DDocas(hit3DList, seedPca);
104 
105  // Use this info to sort the hits along the principle axis - note by absolute value of arc length
106  hit3DList.sort(SeedFinderAlgBase::Sort3DHitsByAbsArcLen3D());
107 
108  // Now translate the seedCenter by the arc len to the first hit
109  double seedDir[3] = {seedPca.getEigenVectors().row(2)(0),
110  seedPca.getEigenVectors().row(2)(1),
111  seedPca.getEigenVectors().row(2)(2)};
112  double seedStart[3] = {
113  hit3DList.front()->getX(), hit3DList.front()->getY(), hit3DList.front()->getZ()};
114 
115  if (hit3DList.size() > 10) {
116  TVector3 newSeedPos;
117  TVector3 newSeedDir;
118  double chiDOF;
119 
120  LineFit2DHits(hit3DList, seedStart[0], newSeedPos, newSeedDir, chiDOF);
121 
122  if (chiDOF > 0.) {
123  // check angles between new/old directions
124  double cosAng = seedDir[0] * newSeedDir[0] + seedDir[1] * newSeedDir[1] +
125  seedDir[2] * newSeedDir[2];
126 
127  if (cosAng < 0.) newSeedDir *= -1.;
128 
129  seedStart[0] = newSeedPos[0];
130  seedStart[1] = newSeedPos[1];
131  seedStart[2] = newSeedPos[2];
132  seedDir[0] = newSeedDir[0];
133  seedDir[1] = newSeedDir[1];
134  seedDir[2] = newSeedDir[2];
135  }
136  }
137 
138  for (const auto& hit3D : hit3DList)
139  hit3D->setStatusBit(0x40000000);
140 
141  seedHitPairVec.emplace_back(std::pair<recob::Seed, reco::HitPairListPtr>(
142  recob::Seed(seedStart, seedDir), hit3DList));
143 
144  foundGoodSeed = true;
145  }
146  }
147  }
148 
149  return foundGoodSeed;
150  }
151 
152  bool
154  reco::PrincipalComponents& seedPca) const
155  {
156  bool foundGoodSeedHits(false);
157 
158  // Use a set to count 2D hits by keeping only a single copy
159  std::set<const reco::ClusterHit2D*> hit2DSet;
160 
161  // Look for numSeed2DHits which are continuous
162  double lastArcLen = hit3DList.front()->getArclenToPoca();
163 
164  reco::HitPairListPtr::iterator startItr = hit3DList.begin();
165  reco::HitPairListPtr::iterator lastItr = hit3DList.begin();
166 
167  while (++lastItr != hit3DList.end()) {
168  const reco::ClusterHit3D* hit3D = *lastItr;
169  double arcLen = hit3D->getArclenToPoca();
170 
171  if (fabs(arcLen - lastArcLen) > m_gapDistance) {
172  startItr = lastItr;
173  hit2DSet.clear();
174  }
175 
176  for (const auto& hit2D : hit3D->getHits()) {
177  hit2DSet.insert(hit2D);
178  }
179 
180  if (hit2DSet.size() > m_numSeed2DHits) break;
181 
182  lastArcLen = arcLen;
183  }
184 
185  if (hit2DSet.size() > m_numSeed2DHits) {
186  if (startItr != hit3DList.begin()) hit3DList.erase(hit3DList.begin(), startItr);
187  if (lastItr != hit3DList.end()) hit3DList.erase(lastItr, hit3DList.end());
188 
189  // On input, the seedPca will contain the original values so we can recover the original axis now
190  Eigen::Vector3f planeVec0(seedPca.getEigenVectors().row(2));
191 
192  m_pcaAlg.PCAAnalysis_3D(hit3DList, seedPca, true);
193 
194  if (seedPca.getSvdOK()) {
195  // Still looking to point "down"
196  if (seedPca.getEigenVectors().row(2)(1) > 0.) seedPca.flipAxis(0);
197 
198  // Check that the seed PCA we have found is consistent with the input PCA
199  Eigen::Vector3f primarySeedAxis(seedPca.getEigenVectors().row(2));
200 
201  double cosAng = primarySeedAxis.dot(planeVec0);
202 
203  // If the proposed seed axis is not relatively aligned with the input PCA then
204  // we should not be using this method to return seeds. Check that here
205  if (cosAng > m_minAllowedCosAng) foundGoodSeedHits = true;
206  }
207  }
208 
209  return foundGoodSeedHits;
210  }
211 
212  //------------------------------------------------------------------------------
213  void
215  double XOrigin,
216  TVector3& Pos,
217  TVector3& Dir,
218  double& ChiDOF) const
219  {
220  // The following is lifted from Bruce Baller to try to get better
221  // initial parameters for a candidate Seed. It is slightly reworked
222  // which is why it is included here instead of used as is.
223  //
224  // Linear fit using X as the independent variable. Hits to be fitted
225  // are passed in the hits vector in a pair form (X, WireID). The
226  // fitted track position at XOrigin is returned in the Pos vector.
227  // The direction cosines are returned in the Dir vector.
228  //
229  // SVD fit adapted from $ROOTSYS/tutorials/matrix/solveLinear.C
230  // Fit equation is w = A(X)v, where w is a vector of hit wires, A is
231  // a matrix to calculate a track projected to a point at X, and v is
232  // a vector (Yo, Zo, dY/dX, dZ/dX).
233  //
234  // Note: The covariance matrix should also be returned
235  // B. Baller August 2014
236 
237  // assume failure
238  ChiDOF = -1;
239 
240  // Make a set of 2D hits so we consider only unique hits
241  std::set<const reco::ClusterHit2D*> hit2DSet;
242 
243  for (const auto& hit3D : hit3DList) {
244  for (const auto& hit2D : hit3D->getHits())
245  hit2DSet.insert(hit2D);
246  }
247 
248  if (hit2DSet.size() < 4) return;
249 
250  const unsigned int nvars = 4;
251  unsigned int npts = 3 * hit3DList.size();
252 
253  TMatrixD A(npts, nvars);
254  // vector holding the Wire number
255  TVectorD w(npts);
256  unsigned short ninpl[3] = {0};
257  unsigned short nok = 0;
258  unsigned short iht(0), cstat, tpc, ipl;
259  double x, cw, sw, off;
260 
261  // Loop over the 2D hits in the above
262  for (const auto& hit : hit2DSet) {
263  geo::WireID wireID = hit->WireID();
264 
265  cstat = wireID.Cryostat;
266  tpc = wireID.TPC;
267  ipl = wireID.Plane;
268 
269  // get the wire plane offset
270  off = m_geometry->WireCoordinate(0, 0, ipl, tpc, cstat);
271 
272  // get the "cosine-like" component
273  cw = m_geometry->WireCoordinate(1, 0, ipl, tpc, cstat) - off;
274 
275  // the "sine-like" component
276  sw = m_geometry->WireCoordinate(0, 1, ipl, tpc, cstat) - off;
277 
278  x = hit->getXPosition() - XOrigin;
279 
280  A[iht][0] = cw;
281  A[iht][1] = sw;
282  A[iht][2] = cw * x;
283  A[iht][3] = sw * x;
284  w[iht] = wireID.Wire - off;
285 
286  ++ninpl[ipl];
287 
288  // need at least two points in a plane
289  if (ninpl[ipl] == 2) ++nok;
290 
291  iht++;
292  }
293 
294  // need at least 2 planes with at least two points
295  if (nok < 2) return;
296 
297  TDecompSVD svd(A);
298  bool ok;
299  TVectorD tVec = svd.Solve(w, ok);
300 
301  ChiDOF = 0;
302 
303  // not enough points to calculate Chisq
304  if (npts <= 4) return;
305 
306  double ypr, zpr, diff;
307 
308  for (const auto& hit : hit2DSet) {
309  geo::WireID wireID = hit->WireID();
310 
311  cstat = wireID.Cryostat;
312  tpc = wireID.TPC;
313  ipl = wireID.Plane;
314  off = m_geometry->WireCoordinate(0, 0, ipl, tpc, cstat);
315  cw = m_geometry->WireCoordinate(1, 0, ipl, tpc, cstat) - off;
316  sw = m_geometry->WireCoordinate(0, 1, ipl, tpc, cstat) - off;
317  x = hit->getXPosition() - XOrigin;
318  ypr = tVec[0] + tVec[2] * x;
319  zpr = tVec[1] + tVec[3] * x;
320  diff = ypr * cw + zpr * sw - (wireID.Wire - off);
321  ChiDOF += diff * diff;
322  }
323 
324  float werr2 = m_geometry->WirePitch() * m_geometry->WirePitch();
325  ChiDOF /= werr2;
326  ChiDOF /= (float)(npts - 4);
327 
328  double norm = sqrt(1 + tVec[2] * tVec[2] + tVec[3] * tVec[3]);
329  Dir[0] = 1 / norm;
330  Dir[1] = tVec[2] / norm;
331  Dir[2] = tVec[3] / norm;
332 
333  Pos[0] = XOrigin;
334  Pos[1] = tVec[0];
335  Pos[2] = tVec[1];
336 
337  } // TrkLineFit()
338 
339  /*
340  //------------------------------------------------------------------------------
341  void PCASeedFinderAlg::LineFit2DHits(const reco::HitPairListPtr& hit3DList,
342  double XOrigin,
343  TVector3& Pos,
344  TVector3& Dir,
345  double& ChiDOF) const
346  {
347  // Linear fit using X as the independent variable. Hits to be fitted
348  // are passed in the hits vector in a pair form (X, WireID). The
349  // fitted track position at XOrigin is returned in the Pos vector.
350  // The direction cosines are returned in the Dir vector.
351  //
352  // SVD fit adapted from $ROOTSYS/tutorials/matrix/solveLinear.C
353  // Fit equation is w = A(X)v, where w is a vector of hit wires, A is
354  // a matrix to calculate a track projected to a point at X, and v is
355  // a vector (Yo, Zo, dY/dX, dZ/dX).
356  //
357  // Note: The covariance matrix should also be returned
358  // B. Baller August 2014
359 
360  // assume failure
361  ChiDOF = -1;
362 
363  if(hit3DList.size() < 4) return;
364 
365  const unsigned int nvars = 4;
366  unsigned int npts = 3 * hit3DList.size();
367 
368  TMatrixD A(npts, nvars);
369  // vector holding the Wire number
370  TVectorD w(npts);
371  unsigned short ninpl[3] = {0};
372  unsigned short nok = 0;
373  unsigned short iht(0), cstat, tpc, ipl;
374  double x, cw, sw, off;
375 
376  for (const auto& hit3D : hit3DList)
377  {
378  // Inner loop over the 2D hits in the above
379  for (const auto& hit : hit3D->getHits())
380  {
381  geo::WireID wireID = hit->getHit().WireID();
382 
383  cstat = wireID.Cryostat;
384  tpc = wireID.TPC;
385  ipl = wireID.Plane;
386 
387  // get the wire plane offset
388  off = m_geometry->WireCoordinate(0, 0, ipl, tpc, cstat);
389 
390  // get the "cosine-like" component
391  cw = m_geometry->WireCoordinate(1, 0, ipl, tpc, cstat) - off;
392 
393  // the "sine-like" component
394  sw = m_geometry->WireCoordinate(0, 1, ipl, tpc, cstat) - off;
395 
396  x = hit->getXPosition() - XOrigin;
397 
398  A[iht][0] = cw;
399  A[iht][1] = sw;
400  A[iht][2] = cw * x;
401  A[iht][3] = sw * x;
402  w[iht] = wireID.Wire - off;
403 
404  ++ninpl[ipl];
405 
406  // need at least two points in a plane
407  if(ninpl[ipl] == 2) ++nok;
408 
409  iht++;
410  }
411  }
412 
413  // need at least 2 planes with at least two points
414  if(nok < 2) return;
415 
416  TDecompSVD svd(A);
417  bool ok;
418  TVectorD tVec = svd.Solve(w, ok);
419 
420  ChiDOF = 0;
421 
422  // not enough points to calculate Chisq
423  if(npts <= 4) return;
424 
425  double ypr, zpr, diff;
426 
427  for(const auto& hit3D : hit3DList)
428  {
429  for (const auto& hit : hit3D->getHits())
430  {
431  geo::WireID wireID = hit->getHit().WireID();
432 
433  cstat = wireID.Cryostat;
434  tpc = wireID.TPC;
435  ipl = wireID.Plane;
436  off = m_geometry->WireCoordinate(0, 0, ipl, tpc, cstat);
437  cw = m_geometry->WireCoordinate(1, 0, ipl, tpc, cstat) - off;
438  sw = m_geometry->WireCoordinate(0, 1, ipl, tpc, cstat) - off;
439  x = hit->getXPosition() - XOrigin;
440  ypr = tVec[0] + tVec[2] * x;
441  zpr = tVec[1] + tVec[3] * x;
442  diff = ypr * cw + zpr * sw - (wireID.Wire - off);
443  ChiDOF += diff * diff;
444  }
445  }
446 
447  float werr2 = m_geometry->WirePitch() * m_geometry->WirePitch();
448  ChiDOF /= werr2;
449  ChiDOF /= (float)(npts - 4);
450 
451  double norm = sqrt(1 + tVec[2] * tVec[2] + tVec[3] * tVec[3]);
452  Dir[0] = 1 / norm;
453  Dir[1] = tVec[2] / norm;
454  Dir[2] = tVec[3] / norm;
455 
456  Pos[0] = XOrigin;
457  Pos[1] = tVec[0];
458  Pos[2] = tVec[1];
459 
460  } // TrkLineFit()
461  */
462 
463 } // 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.
intermediate_table::iterator iterator
void PCAAnalysis_calc3DDocas(const reco::HitPairListPtr &hitPairVector, const reco::PrincipalComponents &pca) const
bool getHitsAtEnd(reco::HitPairListPtr &hit3DList, reco::PrincipalComponents &seedPca) const
Separate function to find hits at the ends of the input 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
void flipAxis(size_t axis)
Definition: Cluster3D.cxx:208
virtual bool findTrackSeeds(reco::HitPairListPtr &hitPairListPtr, reco::PrincipalComponents &inputPCA, SeedHitPairListPairVec &seedHitMap) const
Given the list of hits this will search for candidate Seed objects and return them.
CryostatID_t Cryostat
Index of cryostat.
Definition: geo_types.h:212
WireID_t Wire
Index of the wire within its plane.
Definition: geo_types.h:580
PrincipalComponentsAlg m_pcaAlg
geo::Length_t WirePitch(geo::PlaneID const &planeid) const
Returns the distance between two consecutive wires.
art framework interface to geometry description
const EigenValues & getEigenValues() const
Definition: Cluster3D.h:246
T get(std::string const &key) const
Definition: ParameterSet.h:271
std::list< const reco::ClusterHit3D * > HitPairListPtr
Definition: Cluster3D.h:335
PlaneID_t Plane
Index of the plane within its TPC.
Definition: geo_types.h:493
auto norm(Vector const &v)
Return norm of the specified vector.
Detector simulation of raw signals on wires.
double m_minAllowedCosAng
The minimum cos(ang) between input and seed axes.
PCASeedFinderAlg(fhicl::ParameterSet const &pset)
Constructor.
static unsigned int reverse(QString &chars, unsigned char *level, unsigned int a, unsigned int b)
Definition: qstring.cpp:11649
float getArclenToPoca() const
Definition: Cluster3D.h:170
T copy(T const &v)
#define A
Definition: memgrp.cpp:38
list x
Definition: train.py:276
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
void LineFit2DHits(const reco::HitPairListPtr &hitList, double XOrigin, TVector3 &Pos, TVector3 &Dir, double &ChiDOF) const
const ClusterHit2DVec & getHits() const
Definition: Cluster3D.h:171
const EigenVectors & getEigenVectors() const
Definition: Cluster3D.h:247