LArSoft  v07_13_02
Liquid Argon Software toolkit - http://larsoft.org/
lar_cluster3d::PrincipalComponentsAlg Class Reference

Cluster3D class. More...

#include "PrincipalComponentsAlg.h"

Public Member Functions

 PrincipalComponentsAlg (fhicl::ParameterSet const &pset)
 Constructor. More...
 
virtual ~PrincipalComponentsAlg ()
 Destructor. More...
 
void reconfigure (fhicl::ParameterSet const &pset)
 a handler for the case where the algorithm control parameters are to be reset More...
 
void PCAAnalysis (const reco::HitPairListPtr &hitPairVector, reco::PrincipalComponents &pca, float doca3DScl=3.) const
 Run the Principal Components Analysis. More...
 
void PCAAnalysis_3D (const reco::HitPairListPtr &hitPairList, reco::PrincipalComponents &pca, bool skeletonOnly=false) const
 
void PCAAnalysis_2D (const reco::HitPairListPtr &hitPairVector, reco::PrincipalComponents &pca, bool updateAvePos=false) const
 
void PCAAnalysis_calc3DDocas (const reco::HitPairListPtr &hitPairVector, const reco::PrincipalComponents &pca) const
 
void PCAAnalysis_calc2DDocas (const reco::Hit2DListPtr &hit2DVector, const reco::PrincipalComponents &pca) const
 
int PCAAnalysis_reject2DOutliers (const reco::HitPairListPtr &hitPairVector, reco::PrincipalComponents &pca, float aveHitDoca) const
 
int PCAAnalysis_reject3DOutliers (const reco::HitPairListPtr &hitPairVector, const reco::PrincipalComponents &pca, float aveHitDoca) const
 

Private Member Functions

void getHit2DPocaToAxis (const Eigen::Vector3f &axisPos, const Eigen::Vector3f &axisDir, const reco::ClusterHit2D *hit2D, Eigen::Vector3f &poca, float &arcLenAxis, float &arcLenWire, float &doca)
 This is used to get the poca, doca and arclen along cluster axis to 2D hit. More...
 

Private Attributes

float m_parallel
 means lines are parallel More...
 
geo::Geometrym_geometry
 
const detinfo::DetectorPropertiesm_detector
 

Detailed Description

Cluster3D class.

Definition at line 39 of file PrincipalComponentsAlg.h.

Constructor & Destructor Documentation

lar_cluster3d::PrincipalComponentsAlg::PrincipalComponentsAlg ( fhicl::ParameterSet const &  pset)

Constructor.

Parameters
pset

Definition at line 35 of file PrincipalComponentsAlg.cxx.

References reconfigure().

36 {
37  this->reconfigure(pset);
38 }
void reconfigure(fhicl::ParameterSet const &pset)
a handler for the case where the algorithm control parameters are to be reset
lar_cluster3d::PrincipalComponentsAlg::~PrincipalComponentsAlg ( )
virtual

Destructor.

Definition at line 42 of file PrincipalComponentsAlg.cxx.

43 {
44 }

Member Function Documentation

void lar_cluster3d::PrincipalComponentsAlg::getHit2DPocaToAxis ( const Eigen::Vector3f &  axisPos,
const Eigen::Vector3f &  axisDir,
const reco::ClusterHit2D hit2D,
Eigen::Vector3f &  poca,
float &  arcLenAxis,
float &  arcLenWire,
float &  doca 
)
private

This is used to get the poca, doca and arclen along cluster axis to 2D hit.

Definition at line 55 of file PrincipalComponentsAlg.cxx.

References detinfo::DetectorProperties::ConvertTicksToX(), geo::CryostatID::Cryostat, d, den, geo::WireGeo::Direction(), e, geo::WireGeo::GetCenter(), reco::ClusterHit2D::getHit(), m_detector, m_geometry, m_parallel, recob::Hit::PeakTime(), geo::PlaneID::Plane, geo::TPCID::TPC, recob::Hit::WireID(), and geo::GeometryCore::WireIDToWireGeo().

62 {
63  // Step one is to set up to determine the point of closest approach of this 2D hit to
64  // the cluster's current axis.
65  // Get this wire's geometry object
66  const geo::WireID& hitID = hit2D->getHit().WireID();
67  const geo::WireGeo& wire_geom = m_geometry->WireIDToWireGeo(hitID);
68 
69  // From this, get the parameters of the line for the wire
70  double wirePos[3] = {0.,0.,0.};
71  Eigen::Vector3f wireDir(wire_geom.Direction().X(),wire_geom.Direction().Y(),wire_geom.Direction().Z());
72 
73  wire_geom.GetCenter(wirePos);
74 
75  // Correct the wire position in x to set to correspond to the drift time
76  float hitPeak(hit2D->getHit().PeakTime());
77 
78  wirePos[0] = m_detector->ConvertTicksToX(hitPeak, hitID.Plane, hitID.TPC, hitID.Cryostat);
79 
80  // Get a vector from the wire position to our cluster's current average position
81  Eigen::Vector3f wVec(axisPos(0)-wirePos[0], axisPos(1)-wirePos[1], axisPos(2)-wirePos[2]);
82 
83  // Get the products we need to compute the arc lengths to the distance of closest approach
84  float a(axisDir.dot(axisDir));
85  float b(axisDir.dot(wireDir));
86  float c(wireDir.dot(wireDir));
87  float d(axisDir.dot(wVec));
88  float e(wireDir.dot(wVec));
89 
90  float den(a*c - b*b);
91 
92  // Parallel lines is a special case
93  if (fabs(den) > m_parallel)
94  {
95  arcLenAxis = (b*e - c*d) / den;
96  arcLenWire = (a*e - b*d) / den;
97  }
98  else
99  {
100  mf::LogDebug("Cluster3D") << "** Parallel lines, need a solution here" << std::endl;
101  arcLenAxis = 0.;
102  arcLenWire = 0.;
103  }
104 
105  // Now get the hit position we'll use for the pca analysis
106  poca = Eigen::Vector3f(wirePos[0] + arcLenWire * wireDir(0),
107  wirePos[1] + arcLenWire * wireDir(1),
108  wirePos[2] + arcLenWire * wireDir(2));
109  Eigen::Vector3f axisPocaPos(axisPos[0] + arcLenAxis * axisDir(0),
110  axisPos[1] + arcLenAxis * axisDir(1),
111  axisPos[2] + arcLenAxis * axisDir(2));
112 
113  float deltaX(poca(0) - axisPocaPos(0));
114  float deltaY(poca(1) - axisPocaPos(1));
115  float deltaZ(poca(2) - axisPocaPos(2));
116  float doca2(deltaX*deltaX + deltaY*deltaY + deltaZ*deltaZ);
117 
118  doca = sqrt(doca2);
119 
120  return;
121 }
float m_parallel
means lines are parallel
Geometry description of a TPC wireThe wire is a single straight segment on a wire plane...
Definition: WireGeo.h:61
geo::WireID WireID() const
Initial tdc tick for hit.
Definition: Hit.h:234
Float_t den
Definition: plot.C:37
CryostatID_t Cryostat
Index of cryostat.
Definition: geo_types.h:130
const recob::Hit & getHit() const
Definition: Cluster3D.h:73
Float_t d
Definition: plot.C:237
PlaneID_t Plane
Index of the plane within its TPC.
Definition: geo_types.h:258
Vector Direction() const
Returns the wire direction as a norm-one vector.
Definition: WireGeo.h:377
virtual double ConvertTicksToX(double ticks, int p, int t, int c) const =0
float PeakTime() const
Time of the signal peak, in tick units.
Definition: Hit.h:219
MaybeLogger_< ELseverityLevel::ELsev_success, false > LogDebug
void GetCenter(double *xyz, double localz=0.0) const
Fills the world coordinate of a point on the wire.
Definition: WireGeo.cxx:68
TPCID_t TPC
Index of the TPC within its cryostat.
Definition: geo_types.h:203
Float_t e
Definition: plot.C:34
const detinfo::DetectorProperties * m_detector
WireGeo const & WireIDToWireGeo(geo::WireID const &wireid) const
Returns the specified wire.
void lar_cluster3d::PrincipalComponentsAlg::PCAAnalysis ( const reco::HitPairListPtr hitPairVector,
reco::PrincipalComponents pca,
float  doca3DScl = 3. 
) const

Run the Principal Components Analysis.

Definition at line 151 of file PrincipalComponentsAlg.cxx.

References reco::PrincipalComponents::getAveHitDoca(), reco::PrincipalComponents::getEigenValues(), reco::PrincipalComponents::getEigenVectors(), reco::PrincipalComponents::getSvdOK(), PCAAnalysis_2D(), PCAAnalysis_3D(), PCAAnalysis_reject2DOutliers(), and PCAAnalysis_reject3DOutliers().

152 {
153  // This is the controlling outside function for running
154  // a Principal Components Analysis on the hits in our
155  // input cluster.
156  // There is a bootstrap process to be followed
157  // 1) Get the initial results from all 3D hits associated
158  // with the cluster
159  // 2) Refine the axis by using only the 2D hits on wires
160 
161  // Get the first axis from all 3D hits
162  PCAAnalysis_3D(hitPairVector, pca);
163 
164  // Make sure first pass was good
165  if (!pca.getSvdOK()) return;
166 
167  // First attempt to refine it using only 2D information
168  reco::PrincipalComponents pcaLoop = pca;
169 
170  PCAAnalysis_2D(hitPairVector, pcaLoop);
171 
172  // If valid result then go to next steps
173  if (pcaLoop.getSvdOK())
174  {
175  // Let's check the angle between the original and the updated axis
176  float cosAngle = pcaLoop.getEigenVectors()[0][0] * pca.getEigenVectors()[0][0]
177  + pcaLoop.getEigenVectors()[0][1] * pca.getEigenVectors()[0][1]
178  + pcaLoop.getEigenVectors()[0][2] * pca.getEigenVectors()[0][2];
179 
180  // Set the scale factor for the outlier rejection
181  float sclFctr(3.);
182 
183  // If we had a significant change then let's do some outlier rejection, etc.
184  if (cosAngle < 1.0) // pretty much everyone takes a turn
185  {
186  int maxIterations(3);
187  float maxRange = 3.*sqrt(pcaLoop.getEigenValues()[1]);
188  float aveDoca = pcaLoop.getAveHitDoca(); // was 0.2
189 
190  maxRange = sclFctr * 0.5*(maxRange+aveDoca); // was std::max(maxRange, aveDoca);
191 
192  int numRejHits = PCAAnalysis_reject2DOutliers(hitPairVector, pcaLoop, maxRange);
193  int totalRejects(numRejHits);
194  int maxRejects(0.4*hitPairVector.size());
195 
196  // Try looping to see if we improve things
197  while(maxIterations-- && numRejHits > 0 && totalRejects < maxRejects)
198  {
199  // Run the PCA
200  PCAAnalysis_2D(hitPairVector, pcaLoop, true);
201 
202  maxRange = sclFctr * 0.5*(3.*sqrt(pcaLoop.getEigenValues()[1])+pcaLoop.getAveHitDoca());
203 
204  numRejHits = PCAAnalysis_reject2DOutliers(hitPairVector, pcaLoop, maxRange);
205  }
206 
207  }
208 
209  // Ok at this stage copy the latest results back into the cluster
210  pca = pcaLoop;
211 
212  // Now we make one last pass through the 3D hits to reject outliers there
213  PCAAnalysis_reject3DOutliers(hitPairVector, pca, doca3DScl * pca.getAveHitDoca());
214  }
215  else pca = pcaLoop;
216 
217  return;
218 }
void PCAAnalysis_2D(const reco::HitPairListPtr &hitPairVector, reco::PrincipalComponents &pca, bool updateAvePos=false) const
int PCAAnalysis_reject3DOutliers(const reco::HitPairListPtr &hitPairVector, const reco::PrincipalComponents &pca, float aveHitDoca) const
bool getSvdOK() const
Definition: Cluster3D.h:226
void PCAAnalysis_3D(const reco::HitPairListPtr &hitPairList, reco::PrincipalComponents &pca, bool skeletonOnly=false) const
int PCAAnalysis_reject2DOutliers(const reco::HitPairListPtr &hitPairVector, reco::PrincipalComponents &pca, float aveHitDoca) const
const float * getEigenValues() const
Definition: Cluster3D.h:228
const float getAveHitDoca() const
Definition: Cluster3D.h:231
const EigenVectors & getEigenVectors() const
Definition: Cluster3D.h:229
void lar_cluster3d::PrincipalComponentsAlg::PCAAnalysis_2D ( const reco::HitPairListPtr hitPairVector,
reco::PrincipalComponents pca,
bool  updateAvePos = false 
) const

Definition at line 356 of file PrincipalComponentsAlg.cxx.

References detinfo::DetectorProperties::ConvertTicksToX(), geo::CryostatID::Cryostat, d, den, geo::WireGeo::Direction(), e, reco::PrincipalComponents::getAvePosition(), geo::WireGeo::GetCenter(), reco::PrincipalComponents::getEigenVectors(), art::left(), m_detector, m_geometry, m_parallel, geo::PlaneID::Plane, art::right(), geo::TPCID::TPC, geo::GeometryCore::WireIDToWireGeo(), x, y, and z.

Referenced by PCAAnalysis().

357 {
358  // Once an axis has been found our goal is to refine it by using only the 2D hits
359  // We'll get 3D information for each of these by using the axis as a reference and use
360  // the point of closest approach as the 3D position
361 
362  // Define elements of our covariance matrix
363  float xi2(0.);
364  float xiyi(0.);
365  float xizi(0.);
366  float yi2(0.);
367  float yizi(0.);
368  float zi2(0.);
369 
370  float aveHitDoca(0.);
371  Eigen::Vector3f avePosUpdate(0.,0.,0.);
372  int nHits(0);
373 
374  // Recover existing line parameters for current cluster
375  const reco::PrincipalComponents& inputPca = pca;
376  Eigen::Vector3f avePosition(inputPca.getAvePosition()[0], inputPca.getAvePosition()[1], inputPca.getAvePosition()[2]);
377  Eigen::Vector3f axisDirVec(inputPca.getEigenVectors()[0][0], inputPca.getEigenVectors()[0][1], inputPca.getEigenVectors()[0][2]);
378 
379  // We float loop here so we can use this method for both the first time through
380  // and a second time through where we re-calculate the mean position
381  // So, we need to keep track of the poca which we do with a float vector
382  std::vector<Eigen::Vector3f> hitPosVec;
383 
384  // Outer loop over 3D hits
385  for (const auto& hit3D : hitPairVector)
386  {
387  // Inner loop over 2D hits
388  for (const auto& hit : hit3D->getHits())
389  {
390  // Step one is to set up to determine the point of closest approach of this 2D hit to
391  // the cluster's current axis.
392  // Get this wire's geometry object
393  const geo::WireID& hitID = hit->getHit().WireID();
394  const geo::WireGeo& wire_geom = m_geometry->WireIDToWireGeo(hitID);
395 
396  // From this, get the parameters of the line for the wire
397  double wirePosArr[3] = {0.,0.,0.};
398  wire_geom.GetCenter(wirePosArr);
399 
400  Eigen::Vector3f wireCenter(wirePosArr[0], wirePosArr[1], wirePosArr[2]);
401  Eigen::Vector3f wireDirVec(wire_geom.Direction().X(),wire_geom.Direction().Y(),wire_geom.Direction().Z());
402 
403  // Correct the wire position in x to set to correspond to the drift time
404  float hitPeak(hit->getHit().PeakTime());
405 
406  Eigen::Vector3f wirePos(m_detector->ConvertTicksToX(hitPeak, hitID.Plane, hitID.TPC, hitID.Cryostat), wireCenter[1], wireCenter[2]);
407 
408  // Compute the wire plane normal for this view
409  Eigen::Vector3f xAxis(1.,0.,0.);
410  Eigen::Vector3f planeNormal = xAxis.cross(wireDirVec); // This gives a normal vector in +z for a Y wire
411 
412  float docaInPlane(wirePos[0] - avePosition[0]);
413  float arcLenToPlane(0.);
414  float cosAxisToPlaneNormal = axisDirVec.dot(planeNormal);
415 
416  Eigen::Vector3f axisPlaneIntersection = wirePos;
417  Eigen::Vector3f hitPosTVec = wirePos;
418 
419  if (fabs(cosAxisToPlaneNormal) > 0.)
420  {
421  Eigen::Vector3f deltaPos = wirePos - avePosition;
422 
423  arcLenToPlane = deltaPos.dot(planeNormal) / cosAxisToPlaneNormal;
424  axisPlaneIntersection = avePosition + arcLenToPlane * axisDirVec;
425  docaInPlane = wirePos[0] - axisPlaneIntersection[0];
426 
427  Eigen::Vector3f axisToInter = axisPlaneIntersection - wirePos;
428  float arcLenToDoca = axisToInter.dot(wireDirVec);
429 
430  hitPosTVec += arcLenToDoca * wireDirVec;
431  }
432 
433  // Get a vector from the wire position to our cluster's current average position
434  Eigen::Vector3f wVec = avePosition - wirePos;
435 
436  // Get the products we need to compute the arc lengths to the distance of closest approach
437  float a(axisDirVec.dot(axisDirVec));
438  float b(axisDirVec.dot(wireDirVec));
439  float c(wireDirVec.dot(wireDirVec));
440  float d(axisDirVec.dot(wVec));
441  float e(wireDirVec.dot(wVec));
442 
443  float den(a*c - b*b);
444  float arcLen1(0.);
445  float arcLen2(0.);
446 
447  // Parallel lines is a special case
448  if (fabs(den) > m_parallel)
449  {
450  arcLen1 = (b*e - c*d) / den;
451  arcLen2 = (a*e - b*d) / den;
452  }
453  else
454  {
455  mf::LogDebug("Cluster3D") << "** Parallel lines, need a solution here" << std::endl;
456  break;
457  }
458 
459  // Now get the hit position we'll use for the pca analysis
460  //float hitPos[] = {wirePos[0] + arcLen2 * wireDirVec[0],
461  // wirePos[1] + arcLen2 * wireDirVec[1],
462  // wirePos[2] + arcLen2 * wireDirVec[2]};
463  //float axisPos[] = {avePosition[0] + arcLen1 * axisDirVec[0],
464  // avePosition[1] + arcLen1 * axisDirVec[1],
465  // avePosition[2] + arcLen1 * axisDirVec[2]};
466  Eigen::Vector3f hitPos = wirePos + arcLen2 * wireDirVec;
467  Eigen::Vector3f axisPos = avePosition + arcLen1 * axisDirVec;
468  float deltaX = hitPos(0) - axisPos(0);
469  float deltaY = hitPos(1) - axisPos(1);
470  float deltaZ = hitPos(2) - axisPos(2);
471  float doca2 = deltaX*deltaX + deltaY*deltaY + deltaZ*deltaZ;
472  float doca = sqrt(doca2);
473 
474  docaInPlane = doca;
475 
476  aveHitDoca += fabs(docaInPlane);
477 
478  //Eigen::Vector3f deltaPos = hitPos - hitPosTVec;
479  //float deltaDoca = doca - docaInPlane;
480 
481  //if (fabs(deltaPos[0]) > 1. || fabs(deltaPos[1]) > 1. || fabs(deltaPos[2]) > 1. || fabs(deltaDoca) > 2.)
482  //{
483  // std::cout << "**************************************************************************************" << std::endl;
484  // std::cout << "Diff in doca: " << deltaPos[0] << "," << deltaPos[1] << "," << deltaPos[2] << ", deltaDoca: " << deltaDoca << std::endl;
485  // std::cout << "-- HitPosTVec: " << hitPosTVec[0] << "," << hitPosTVec[1] << "," << hitPosTVec[2] << std::endl;
486  // std::cout << "-- hitPos: " << hitPos[0] << "," << hitPos[1] << "," << hitPos[2] << std::endl;
487  // std::cout << "-- WirePos: " << wirePos[0] << "," << wirePos[1] << "," << wirePos[2] << std::endl;
488  //}
489 
490  // Set the hit's doca and arclen
491  hit->setDocaToAxis(fabs(docaInPlane));
492  hit->setArcLenToPoca(arcLenToPlane);
493 
494  // If this point is considered an outlier then we skip
495  // the accumulation for average position and covariance
496  if (hit->getStatusBits() & 0x80)
497  {
498  continue;
499  }
500 
501  //hitPosTVec = hitPos;
502 
503  avePosUpdate += hitPosTVec;
504 
505  hitPosVec.push_back(hitPosTVec);
506 
507  nHits++;
508  }
509  }
510 
511  // Get updated average position
512  avePosUpdate[0] /= float(nHits);
513  avePosUpdate[1] /= float(nHits);
514  avePosUpdate[2] /= float(nHits);
515 
516  // Get the average hit doca
517  aveHitDoca /= float(nHits);
518 
519  if (updateAvePos)
520  {
521  avePosition = avePosUpdate;
522  }
523 
524  // Now loop through the hits and build out the covariance matrix
525  for(auto& hitPos : hitPosVec)
526  {
527  // And increment the values in the covariance matrix
528  float x = hitPos[0] - avePosition[0];
529  float y = hitPos[1] - avePosition[1];
530  float z = hitPos[2] - avePosition[2];
531 
532  xi2 += x * x;
533  xiyi += x * y;
534  xizi += x * z;
535  yi2 += y * y;
536  yizi += y * z;
537  zi2 += z * z;
538  }
539 
540  // Accumulation done, now do the actual work
541 
542  // Using Eigen package
543  Eigen::Matrix3f sig;
544 
545  sig << xi2, xiyi, xizi,
546  xiyi, yi2, yizi,
547  xizi, yizi, zi2;
548 
549  sig *= 1./(nHits - 1);
550 
551  Eigen::SelfAdjointEigenSolver<Eigen::Matrix3f> eigenMat(sig);
552 
553  if (eigenMat.info() == Eigen::ComputationInfo::Success)
554  {
555  using eigenValColPair = std::pair<float,size_t>;
556  std::vector<eigenValColPair> eigenValColVec;
557 
558  eigenValColVec.push_back(eigenValColPair(eigenMat.eigenvalues()(0),0));
559  eigenValColVec.push_back(eigenValColPair(eigenMat.eigenvalues()(1),1));
560  eigenValColVec.push_back(eigenValColPair(eigenMat.eigenvalues()(2),2));
561 
562  std::sort(eigenValColVec.begin(),eigenValColVec.end(),[](const eigenValColPair& left, const eigenValColPair& right){return left.first > right.first;});
563 
564  // Now copy output
565  // Get the eigen values
566  float recobEigenVals[] = {eigenValColVec[0].first, eigenValColVec[1].first, eigenValColVec[2].first};
567 
568  // Grab the principle axes
570  Eigen::Matrix3f eigenVecs(eigenMat.eigenvectors());
571 
572  for(const auto& pair : eigenValColVec)
573  {
574  std::vector<float> tempVec = {eigenVecs(0,pair.second),eigenVecs(1,pair.second),eigenVecs(2,pair.second)};
575  recobEigenVecs.push_back(tempVec);
576  }
577 
578  // Save the average position
579  float avePosToSave[] = {float(avePosition[0]),float(avePosition[1]),float(avePosition[2])};
580 
581  // Store away
582  pca = reco::PrincipalComponents(true, nHits, recobEigenVals, recobEigenVecs, avePosToSave, aveHitDoca);
583  }
584  else
585  {
586  mf::LogDebug("Cluster3D") << "PCA decompose failure, numPairs = " << nHits << std::endl;
588  }
589 
590  return;
591 }
Float_t x
Definition: compare.C:6
float m_parallel
means lines are parallel
Geometry description of a TPC wireThe wire is a single straight segment on a wire plane...
Definition: WireGeo.h:61
constexpr auto const & right(const_AssnsIter< L, R, D, Dir > const &a, const_AssnsIter< L, R, D, Dir > const &b)
Definition: AssnsIter.h:112
Float_t den
Definition: plot.C:37
Float_t y
Definition: compare.C:6
Double_t z
Definition: plot.C:279
CryostatID_t Cryostat
Index of cryostat.
Definition: geo_types.h:130
const float * getAvePosition() const
Definition: Cluster3D.h:230
Float_t d
Definition: plot.C:237
std::vector< std::vector< float > > EigenVectors
Definition: Cluster3D.h:209
PlaneID_t Plane
Index of the plane within its TPC.
Definition: geo_types.h:258
Vector Direction() const
Returns the wire direction as a norm-one vector.
Definition: WireGeo.h:377
Detector simulation of raw signals on wires.
virtual double ConvertTicksToX(double ticks, int p, int t, int c) const =0
constexpr auto const & left(const_AssnsIter< L, R, D, Dir > const &a, const_AssnsIter< L, R, D, Dir > const &b)
Definition: AssnsIter.h:104
MaybeLogger_< ELseverityLevel::ELsev_success, false > LogDebug
void GetCenter(double *xyz, double localz=0.0) const
Fills the world coordinate of a point on the wire.
Definition: WireGeo.cxx:68
TPCID_t TPC
Index of the TPC within its cryostat.
Definition: geo_types.h:203
Float_t e
Definition: plot.C:34
const detinfo::DetectorProperties * m_detector
const EigenVectors & getEigenVectors() const
Definition: Cluster3D.h:229
WireGeo const & WireIDToWireGeo(geo::WireID const &wireid) const
Returns the specified wire.
void lar_cluster3d::PrincipalComponentsAlg::PCAAnalysis_3D ( const reco::HitPairListPtr hitPairList,
reco::PrincipalComponents pca,
bool  skeletonOnly = false 
) const

Definition at line 220 of file PrincipalComponentsAlg.cxx.

References art::left(), max, art::right(), reco::ClusterHit3D::SKELETONHIT, weight, x, y, and z.

Referenced by lar_cluster3d::ClusterPathFinder::breakIntoTinyBits(), lar_cluster3d::VoronoiPathFinder::breakIntoTinyBits(), lar_cluster3d::HoughSeedFinderAlg::buildSeed(), lar_cluster3d::MinSpanTreeAlg::CheckHitSorting(), lar_cluster3d::ConvexHullPathFinder::completeCandidateCluster(), lar_cluster3d::ClusterParamsBuilder::FillClusterParams(), lar_cluster3d::MinSpanTreeAlg::FindBestPathInCluster(), lar_cluster3d::HoughSeedFinderAlg::findHitGaps(), lar_cluster3d::HoughSeedFinderAlg::findTrackSeeds(), lar_cluster3d::ParallelHitsSeedFinderAlg::findTrackSeeds(), lar_cluster3d::PCASeedFinderAlg::getHitsAtEnd(), lar_cluster3d::VoronoiPathFinder::makeCandidateCluster(), lar_cluster3d::ClusterMergeAlg::mergeClusters(), PCAAnalysis(), and lar_cluster3d::Cluster3D::splitClustersWithHough().

221 {
222  // We want to run a PCA on the input TkrVecPoints...
223  // The steps are:
224  // 1) do a mean normalization of the input vec points
225  // 2) compute the covariance matrix
226  // 3) run the SVD
227  // 4) extract the eigen vectors and values
228  // see what happens
229 
230  // Run through the HitPairList and get the mean position of all the hits
231  float meanPos[] = {0.,0.,0.};
232  float meanWeightSum(0.);
233  int numPairsInt(0);
234 
235 // const float minimumDeltaPeakSig(0.00001);
236  float minimumDeltaPeakSig(0.00001);
237 
238  // Want to use the hit "chi square" to weight the hits but we need to put a lower limit on its value
239  // to prevent a few hits being over counted.
240  // This is a bit experimental until we can evaluate the cost (time to calculate) vs the benefit
241  // (better fits)..
242  std::vector<float> hitChiSquareVec;
243 
244  hitChiSquareVec.resize(hitPairVector.size());
245 
246  std::transform(hitPairVector.begin(),hitPairVector.end(),hitChiSquareVec.begin(),[](const auto& hit){return hit->getHitChiSquare();});
247  std::sort(hitChiSquareVec.begin(),hitChiSquareVec.end());
248 
249  size_t numToKeep = 0.8 * hitChiSquareVec.size();
250 
251  hitChiSquareVec.resize(numToKeep);
252 
253  float aveValue = std::accumulate(hitChiSquareVec.begin(),hitChiSquareVec.end(),float(0.)) / float(hitChiSquareVec.size());
254  float rms = std::sqrt(std::inner_product(hitChiSquareVec.begin(),hitChiSquareVec.end(), hitChiSquareVec.begin(), 0.,std::plus<>(),[aveValue](const auto& left,const auto& right){return (left - aveValue) * (right - aveValue);}) / float(hitChiSquareVec.size()));
255 
256  minimumDeltaPeakSig = std::max(minimumDeltaPeakSig, aveValue - rms);
257 
258 // std::cout << "===>> Calculating PCA, ave chiSquare: " << aveValue << ", rms: " << rms << ", cut: " << minimumDeltaPeakSig << std::endl;
259 
260  for (const auto& hit : hitPairVector)
261  {
262  if (skeletonOnly && !((hit->getStatusBits() & reco::ClusterHit3D::SKELETONHIT) == reco::ClusterHit3D::SKELETONHIT)) continue;
263 
264  // Weight the hit by the peak time difference significance
265  float weight = std::max(minimumDeltaPeakSig, hit->getHitChiSquare()); //hit->getDeltaPeakTime()); ///hit->getSigmaPeakTime());
266 
267  meanPos[0] += hit->getPosition()[0] * weight;
268  meanPos[1] += hit->getPosition()[1] * weight;
269  meanPos[2] += hit->getPosition()[2] * weight;
270  numPairsInt++;
271 
272  meanWeightSum += weight;
273  }
274 
275  meanPos[0] /= meanWeightSum;
276  meanPos[1] /= meanWeightSum;
277  meanPos[2] /= meanWeightSum;
278 
279  // Define elements of our covariance matrix
280  float xi2(0.);
281  float xiyi(0.);
282  float xizi(0.0);
283  float yi2(0.0);
284  float yizi(0.0);
285  float zi2(0.);
286  float weightSum(0.);
287 
288  // Back through the hits to build the matrix
289  for (const auto& hit : hitPairVector)
290  {
291  if (skeletonOnly && !((hit->getStatusBits() & reco::ClusterHit3D::SKELETONHIT) == reco::ClusterHit3D::SKELETONHIT)) continue;
292 
293  float weight = 1. / std::max(minimumDeltaPeakSig, hit->getHitChiSquare()); //hit->getDeltaPeakTime()); ///hit->getSigmaPeakTime());
294  float x = (hit->getPosition()[0] - meanPos[0]) * weight;
295  float y = (hit->getPosition()[1] - meanPos[1]) * weight;
296  float z = (hit->getPosition()[2] - meanPos[2]) * weight;
297 
298  weightSum += weight*weight;
299 
300  xi2 += x * x;
301  xiyi += x * y;
302  xizi += x * z;
303  yi2 += y * y;
304  yizi += y * z;
305  zi2 += z * z;
306  }
307 
308  // Using Eigen package
309  Eigen::Matrix3f sig;
310 
311  sig << xi2, xiyi, xizi,
312  xiyi, yi2, yizi,
313  xizi, yizi, zi2;
314 
315  sig *= 1./weightSum;
316 
317  Eigen::SelfAdjointEigenSolver<Eigen::Matrix3f> eigenMat(sig);
318 
319  if (eigenMat.info() == Eigen::ComputationInfo::Success)
320  {
321  using eigenValColPair = std::pair<float,size_t>;
322  std::vector<eigenValColPair> eigenValColVec;
323 
324  eigenValColVec.push_back(eigenValColPair(eigenMat.eigenvalues()(0),0));
325  eigenValColVec.push_back(eigenValColPair(eigenMat.eigenvalues()(1),1));
326  eigenValColVec.push_back(eigenValColPair(eigenMat.eigenvalues()(2),2));
327 
328  std::sort(eigenValColVec.begin(),eigenValColVec.end(),[](const eigenValColPair& left, const eigenValColPair& right){return left.first > right.first;});
329 
330  // Now copy output
331  // Get the eigen values
332  float recobEigenVals[] = {eigenValColVec[0].first, eigenValColVec[1].first, eigenValColVec[2].first};
333 
334  // Grab the principle axes
336  Eigen::Matrix3f eigenVecs(eigenMat.eigenvectors());
337 
338  for(const auto& pair : eigenValColVec)
339  {
340  std::vector<float> tempVec = {eigenVecs(0,pair.second),eigenVecs(1,pair.second),eigenVecs(2,pair.second)};
341  recobEigenVecs.push_back(tempVec);
342  }
343 
344  // Store away
345  pca = reco::PrincipalComponents(true, numPairsInt, recobEigenVals, recobEigenVecs, meanPos);
346  }
347  else
348  {
349  mf::LogDebug("Cluster3D") << "PCA decompose failure, numPairs = " << numPairsInt << std::endl;
351  }
352 
353  return;
354 }
Float_t x
Definition: compare.C:6
constexpr auto const & right(const_AssnsIter< L, R, D, Dir > const &a, const_AssnsIter< L, R, D, Dir > const &b)
Definition: AssnsIter.h:112
Float_t y
Definition: compare.C:6
Double_t z
Definition: plot.C:279
Int_t max
Definition: plot.C:27
std::vector< std::vector< float > > EigenVectors
Definition: Cluster3D.h:209
Detector simulation of raw signals on wires.
constexpr auto const & left(const_AssnsIter< L, R, D, Dir > const &a, const_AssnsIter< L, R, D, Dir > const &b)
Definition: AssnsIter.h:104
double weight
Definition: plottest35.C:25
MaybeLogger_< ELseverityLevel::ELsev_success, false > LogDebug
Hit is a "skeleton" hit.
Definition: Cluster3D.h:93
void lar_cluster3d::PrincipalComponentsAlg::PCAAnalysis_calc2DDocas ( const reco::Hit2DListPtr hit2DVector,
const reco::PrincipalComponents pca 
) const

Definition at line 645 of file PrincipalComponentsAlg.cxx.

References geo::WireGeo::Direction(), reco::PrincipalComponents::getAvePosition(), geo::WireGeo::GetCenter(), reco::PrincipalComponents::getEigenValues(), reco::PrincipalComponents::getEigenVectors(), geo::WireGeo::HalfL(), m_geometry, min, reco::PrincipalComponents::setAveHitDoca(), and geo::GeometryCore::WireIDToWireGeo().

647 {
648  // Our mission, should we choose to accept it, is to scan through the 2D hits and reject
649  // any outliers. Basically, any hit outside a scaled range of the average doca from the
650  // first pass is marked by setting the bit in the status word.
651 
652  // We'll need the current PCA axis to determine doca and arclen
653  Eigen::Vector3f avePosition(pca.getAvePosition()[0], pca.getAvePosition()[1], pca.getAvePosition()[2]);
654  Eigen::Vector3f axisDirVec(pca.getEigenVectors()[0][0], pca.getEigenVectors()[0][1], pca.getEigenVectors()[0][2]);
655 
656  // Recover the principle eigen value for range constraints
657  float maxArcLen = 4.*sqrt(pca.getEigenValues()[0]);
658 
659  // We want to keep track of the average
660  float aveHitDoca(0.);
661 
662  // Outer loop over views
663  for (const auto* hit : hit2DListPtr)
664  {
665  // Step one is to set up to determine the point of closest approach of this 2D hit to
666  // the cluster's current axis. We do that by finding the point of intersection of the
667  // cluster's axis with a plane defined by the wire the hit is associated with.
668  // Get this wire's geometry object
669  const geo::WireID& hitID = hit->getHit().WireID();
670  const geo::WireGeo& wire_geom = m_geometry->WireIDToWireGeo(hitID);
671 
672  // From this, get the parameters of the line for the wire
673  double wirePosArr[3] = {0.,0.,0.};
674  wire_geom.GetCenter(wirePosArr);
675 
676  Eigen::Vector3f wireCenter(wirePosArr[0], wirePosArr[1], wirePosArr[2]);
677  Eigen::Vector3f wireDirVec(wire_geom.Direction().X(),wire_geom.Direction().Y(),wire_geom.Direction().Z());
678 
679  // Correct the wire position in x to set to correspond to the drift time
680  Eigen::Vector3f wirePos(hit->getXPosition(), wireCenter[1], wireCenter[2]);
681 
682  // Compute the wire plane normal for this view
683  Eigen::Vector3f xAxis(1.,0.,0.);
684  Eigen::Vector3f planeNormal = xAxis.cross(wireDirVec); // This gives a normal vector in +z for a Y wire
685 
686  float arcLenToPlane(0.);
687  float docaInPlane(wirePos[0] - avePosition[0]);
688  float cosAxisToPlaneNormal = axisDirVec.dot(planeNormal);
689 
690  Eigen::Vector3f axisPlaneIntersection = wirePos;
691 
692  // If current cluster axis is not parallel to wire plane then find intersection point
693  if (fabs(cosAxisToPlaneNormal) > 0.)
694  {
695  Eigen::Vector3f deltaPos = wirePos - avePosition;
696 
697  arcLenToPlane = std::min(float(deltaPos.dot(planeNormal) / cosAxisToPlaneNormal), maxArcLen);
698  axisPlaneIntersection = avePosition + arcLenToPlane * axisDirVec;
699 
700  Eigen::Vector3f axisToInter = axisPlaneIntersection - wirePos;
701  float arcLenToDoca = axisToInter.dot(wireDirVec);
702 
703  // If the arc length along the wire to the poca is outside the TPC then reset
704  if (fabs(arcLenToDoca) > wire_geom.HalfL()) arcLenToDoca = wire_geom.HalfL();
705 
706  // If we were successful in getting to the wire plane then the doca is simply the
707  // difference in x coordinates... but we hvae to worry about the special cases so
708  // we calculate a 3D doca based on arclengths above...
709  Eigen::Vector3f docaVec = axisPlaneIntersection - (wirePos + arcLenToDoca * wireDirVec);
710  docaInPlane = docaVec.norm();
711  }
712 
713  aveHitDoca += fabs(docaInPlane);
714 
715  // Set the hit's doca and arclen
716  hit->setDocaToAxis(fabs(docaInPlane));
717  hit->setArcLenToPoca(arcLenToPlane);
718  }
719 
720  // Compute the average and store
721  aveHitDoca /= float(hit2DListPtr.size());
722 
723  pca.setAveHitDoca(aveHitDoca);
724 
725  return;
726 }
Geometry description of a TPC wireThe wire is a single straight segment on a wire plane...
Definition: WireGeo.h:61
void setAveHitDoca(double doca) const
Definition: Cluster3D.h:234
const float * getAvePosition() const
Definition: Cluster3D.h:230
Vector Direction() const
Returns the wire direction as a norm-one vector.
Definition: WireGeo.h:377
Detector simulation of raw signals on wires.
double HalfL() const
Returns half the length of the wire [cm].
Definition: WireGeo.h:106
const float * getEigenValues() const
Definition: Cluster3D.h:228
Int_t min
Definition: plot.C:26
void GetCenter(double *xyz, double localz=0.0) const
Fills the world coordinate of a point on the wire.
Definition: WireGeo.cxx:68
const EigenVectors & getEigenVectors() const
Definition: Cluster3D.h:229
WireGeo const & WireIDToWireGeo(geo::WireID const &wireid) const
Returns the specified wire.
void lar_cluster3d::PrincipalComponentsAlg::PCAAnalysis_calc3DDocas ( const reco::HitPairListPtr hitPairVector,
const reco::PrincipalComponents pca 
) const

Definition at line 593 of file PrincipalComponentsAlg.cxx.

References reco::PrincipalComponents::getAvePosition(), reco::PrincipalComponents::getEigenVectors(), and reco::PrincipalComponents::setAveHitDoca().

Referenced by lar_cluster3d::ClusterPathFinder::breakIntoTinyBits(), lar_cluster3d::VoronoiPathFinder::breakIntoTinyBits(), lar_cluster3d::HoughSeedFinderAlg::buildSeed(), lar_cluster3d::HoughSeedFinderAlg::findHitGaps(), lar_cluster3d::PCASeedFinderAlg::findTrackSeeds(), lar_cluster3d::ParallelHitsSeedFinderAlg::findTrackSeeds(), lar_cluster3d::Cluster3D::splitClustersWithHough(), lar_cluster3d::ConvexHullPathFinder::subDivideCluster(), and lar_cluster3d::VoronoiPathFinder::subDivideCluster().

595 {
596  // Our mission, should we choose to accept it, is to scan through the 2D hits and reject
597  // any outliers. Basically, any hit outside a scaled range of the average doca from the
598  // first pass is marked by setting the bit in the status word.
599 
600  // We'll need the current PCA axis to determine doca and arclen
601  Eigen::Vector3f avePosition(pca.getAvePosition()[0], pca.getAvePosition()[1], pca.getAvePosition()[2]);
602  Eigen::Vector3f axisDirVec(pca.getEigenVectors()[0][0], pca.getEigenVectors()[0][1], pca.getEigenVectors()[0][2]);
603 
604  // We want to keep track of the average
605  float aveDoca3D(0.);
606 
607  // Outer loop over views
608  for (const auto* clusterHit3D : hitPairVector)
609  {
610  // Always reset the existing status bit
611  clusterHit3D->clearStatusBits(0x80);
612 
613  // Now we need to calculate the doca and poca...
614  // Start by getting this hits position
615  Eigen::Vector3f clusPos(clusterHit3D->getPosition()[0],clusterHit3D->getPosition()[1],clusterHit3D->getPosition()[2]);
616 
617  // Form a TVector from this to the cluster average position
618  Eigen::Vector3f clusToHitVec = clusPos - avePosition;
619 
620  // With this we can get the arclength to the doca point
621  float arclenToPoca = clusToHitVec.dot(axisDirVec);
622 
623  // Get the coordinates along the axis for this point
624  Eigen::Vector3f docaPos = avePosition + arclenToPoca * axisDirVec;
625 
626  // Now get doca and poca
627  Eigen::Vector3f docaPosToClusPos = clusPos - docaPos;
628  float docaToAxis = docaPosToClusPos.norm();
629 
630  aveDoca3D += docaToAxis;
631 
632  // Ok, set the values in the hit
633  clusterHit3D->setDocaToAxis(docaToAxis);
634  clusterHit3D->setArclenToPoca(arclenToPoca);
635  }
636 
637  // Compute the average and store
638  aveDoca3D /= float(hitPairVector.size());
639 
640  pca.setAveHitDoca(aveDoca3D);
641 
642  return;
643 }
void setAveHitDoca(double doca) const
Definition: Cluster3D.h:234
const float * getAvePosition() const
Definition: Cluster3D.h:230
const EigenVectors & getEigenVectors() const
Definition: Cluster3D.h:229
int lar_cluster3d::PrincipalComponentsAlg::PCAAnalysis_reject2DOutliers ( const reco::HitPairListPtr hitPairVector,
reco::PrincipalComponents pca,
float  aveHitDoca 
) const

Definition at line 728 of file PrincipalComponentsAlg.cxx.

Referenced by PCAAnalysis().

731 {
732  // Our mission, should we choose to accept it, is to scan through the 2D hits and reject
733  // any outliers. Basically, any hit outside a scaled range of the average doca from the
734  // first pass is marked by setting the bit in the status word.
735 
736  // First get the average doca scaled by some appropriate factor
737  int numRejHits(0);
738 
739  // Outer loop over views
740  for (const auto& hit3D : hitPairVector)
741  {
742  // Inner loop over hits in this view
743  for (const auto& hit : hit3D->getHits())
744  {
745  // Always reset the existing status bit
746  hit->clearStatusBits(0x80);
747 
748  if (hit->getDocaToAxis() > maxDocaAllowed)
749  {
750  hit->setStatusBit(0x80);
751  numRejHits++;
752  }
753  }
754  }
755 
756  return numRejHits;
757 }
Detector simulation of raw signals on wires.
int lar_cluster3d::PrincipalComponentsAlg::PCAAnalysis_reject3DOutliers ( const reco::HitPairListPtr hitPairVector,
const reco::PrincipalComponents pca,
float  aveHitDoca 
) const

Definition at line 759 of file PrincipalComponentsAlg.cxx.

References reco::PrincipalComponents::getAvePosition(), and reco::PrincipalComponents::getEigenVectors().

Referenced by PCAAnalysis().

762 {
763  // Our mission, should we choose to accept it, is to scan through the 2D hits and reject
764  // any outliers. Basically, any hit outside a scaled range of the average doca from the
765  // first pass is marked by setting the bit in the status word.
766 
767  // First get the average doca scaled by some appropriate factor
768  int numRejHits(0);
769 
770  // We'll need the current PCA axis to determine doca and arclen
771  Eigen::Vector3f avePosition(pca.getAvePosition()[0], pca.getAvePosition()[1], pca.getAvePosition()[2]);
772  Eigen::Vector3f axisDirVec(pca.getEigenVectors()[0][0], pca.getEigenVectors()[0][1], pca.getEigenVectors()[0][2]);
773 
774  // Outer loop over views
775  for (const auto* clusterHit3D : hitPairVector)
776  {
777  // Always reset the existing status bit
778  clusterHit3D->clearStatusBits(0x80);
779 
780  // Now we need to calculate the doca and poca...
781  // Start by getting this hits position
782  Eigen::Vector3f clusPos(clusterHit3D->getPosition()[0],clusterHit3D->getPosition()[1],clusterHit3D->getPosition()[2]);
783 
784  // Form a TVector from this to the cluster average position
785  Eigen::Vector3f clusToHitVec = clusPos - avePosition;
786 
787  // With this we can get the arclength to the doca point
788  float arclenToPoca = clusToHitVec.dot(axisDirVec);
789 
790  // Get the coordinates along the axis for this point
791  Eigen::Vector3f docaPos = avePosition + arclenToPoca * axisDirVec;
792 
793  // Now get doca and poca
794  Eigen::Vector3f docaPosToClusPos = clusPos - docaPos;
795  float docaToAxis = docaPosToClusPos.norm();
796 
797  // Ok, set the values in the hit
798  clusterHit3D->setDocaToAxis(docaToAxis);
799  clusterHit3D->setArclenToPoca(arclenToPoca);
800 
801  // Check to see if this is a keeper
802  if (clusterHit3D->getDocaToAxis() > maxDocaAllowed)
803  {
804  clusterHit3D->setStatusBit(0x80);
805  numRejHits++;
806  }
807  }
808 
809  return numRejHits;
810 }
const float * getAvePosition() const
Definition: Cluster3D.h:230
const EigenVectors & getEigenVectors() const
Definition: Cluster3D.h:229
void lar_cluster3d::PrincipalComponentsAlg::reconfigure ( fhicl::ParameterSet const &  pset)

a handler for the case where the algorithm control parameters are to be reset

Definition at line 46 of file PrincipalComponentsAlg.cxx.

References fhicl::ParameterSet::get(), m_detector, m_geometry, and m_parallel.

Referenced by PrincipalComponentsAlg(), lar_cluster3d::HoughSeedFinderAlg::reconfigure(), lar_cluster3d::PCASeedFinderAlg::reconfigure(), lar_cluster3d::ParallelHitsSeedFinderAlg::reconfigure(), cosmic::CosmicPCAxisTagger::reconfigure(), and lar_cluster3d::Cluster3D::reconfigure().

47 {
49 
50  m_parallel = pset.get<float>("ParallelLines", 0.00001);
51  m_geometry = &*geometry;
52  m_detector = lar::providerFrom<detinfo::DetectorPropertiesService>();
53 }
float m_parallel
means lines are parallel
const detinfo::DetectorProperties * m_detector

Member Data Documentation

const detinfo::DetectorProperties* lar_cluster3d::PrincipalComponentsAlg::m_detector
private

Definition at line 93 of file PrincipalComponentsAlg.h.

Referenced by getHit2DPocaToAxis(), PCAAnalysis_2D(), and reconfigure().

geo::Geometry* lar_cluster3d::PrincipalComponentsAlg::m_geometry
private
float lar_cluster3d::PrincipalComponentsAlg::m_parallel
private

means lines are parallel

Definition at line 90 of file PrincipalComponentsAlg.h.

Referenced by getHit2DPocaToAxis(), PCAAnalysis_2D(), and reconfigure().


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