All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Private Attributes | List of all members
lar_cluster3d::PrincipalComponentsAlg Class Reference

Cluster3D class. More...

#include <PrincipalComponentsAlg.h>

Public Member Functions

 PrincipalComponentsAlg (fhicl::ParameterSet const &pset)
 Constructor. More...
 
void PCAAnalysis (const detinfo::DetectorPropertiesData &detProp, 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 detinfo::DetectorPropertiesData &detProp, 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 Attributes

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

Detailed Description

Cluster3D class.

Definition at line 33 of file PrincipalComponentsAlg.h.

Constructor & Destructor Documentation

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

Constructor.

Parameters
pset

Definition at line 40 of file PrincipalComponentsAlg.cxx.

41  {
42  m_parallel = pset.get<float>("ParallelLines", 0.00001);
43  m_geometry = art::ServiceHandle<geo::Geometry const>{}.get();
44  }
float m_parallel
means lines are parallel

Member Function Documentation

void lar_cluster3d::PrincipalComponentsAlg::PCAAnalysis ( const detinfo::DetectorPropertiesData detProp,
const reco::HitPairListPtr hitPairVector,
reco::PrincipalComponents pca,
float  doca3DScl = 3. 
) const

Run the Principal Components Analysis.

Definition at line 71 of file PrincipalComponentsAlg.cxx.

75  {
76  // This is the controlling outside function for running
77  // a Principal Components Analysis on the hits in our
78  // input cluster.
79  // There is a bootstrap process to be followed
80  // 1) Get the initial results from all 3D hits associated
81  // with the cluster
82  // 2) Refine the axis by using only the 2D hits on wires
83 
84  // Get the first axis from all 3D hits
85  PCAAnalysis_3D(hitPairVector, pca);
86 
87  // Make sure first pass was good
88  if (!pca.getSvdOK()) return;
89 
90  // First attempt to refine it using only 2D information
91  reco::PrincipalComponents pcaLoop = pca;
92 
93  PCAAnalysis_2D(detProp, hitPairVector, pcaLoop);
94 
95  // If valid result then go to next steps
96  if (pcaLoop.getSvdOK()) {
97  // Let's check the angle between the original and the updated axis
98  float cosAngle = pcaLoop.getEigenVectors().row(0) * pca.getEigenVectors().row(0).transpose();
99 
100  // Set the scale factor for the outlier rejection
101  float sclFctr(3.);
102 
103  // If we had a significant change then let's do some outlier rejection, etc.
104  if (cosAngle < 1.0) // pretty much everyone takes a turn
105  {
106  int maxIterations(3);
107  float maxRange = 3. * sqrt(pcaLoop.getEigenValues()[1]);
108  float aveDoca = pcaLoop.getAveHitDoca(); // was 0.2
109 
110  maxRange = sclFctr * 0.5 * (maxRange + aveDoca); // was std::max(maxRange, aveDoca);
111 
112  int numRejHits = PCAAnalysis_reject2DOutliers(hitPairVector, pcaLoop, maxRange);
113  int totalRejects(numRejHits);
114  int maxRejects(0.4 * hitPairVector.size());
115 
116  // Try looping to see if we improve things
117  while (maxIterations-- && numRejHits > 0 && totalRejects < maxRejects) {
118  // Run the PCA
119  PCAAnalysis_2D(detProp, hitPairVector, pcaLoop, true);
120 
121  maxRange =
122  sclFctr * 0.5 * (3. * sqrt(pcaLoop.getEigenValues()[1]) + pcaLoop.getAveHitDoca());
123 
124  numRejHits = PCAAnalysis_reject2DOutliers(hitPairVector, pcaLoop, maxRange);
125  }
126  }
127 
128  // Ok at this stage copy the latest results back into the cluster
129  pca = pcaLoop;
130 
131  // Now we make one last pass through the 3D hits to reject outliers there
132  PCAAnalysis_reject3DOutliers(hitPairVector, pca, doca3DScl * pca.getAveHitDoca());
133  }
134  else
135  pca = pcaLoop;
136 
137  return;
138  }
int PCAAnalysis_reject3DOutliers(const reco::HitPairListPtr &hitPairVector, const reco::PrincipalComponents &pca, float aveHitDoca) const
bool getSvdOK() const
Definition: Cluster3D.h:243
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 EigenValues & getEigenValues() const
Definition: Cluster3D.h:245
void PCAAnalysis_2D(const detinfo::DetectorPropertiesData &detProp, const reco::HitPairListPtr &hitPairVector, reco::PrincipalComponents &pca, bool updateAvePos=false) const
const float getAveHitDoca() const
Definition: Cluster3D.h:248
const EigenVectors & getEigenVectors() const
Definition: Cluster3D.h:246
void lar_cluster3d::PrincipalComponentsAlg::PCAAnalysis_2D ( const detinfo::DetectorPropertiesData detProp,
const reco::HitPairListPtr hitPairVector,
reco::PrincipalComponents pca,
bool  updateAvePos = false 
) const

Definition at line 287 of file PrincipalComponentsAlg.cxx.

291  {
292  // Once an axis has been found our goal is to refine it by using only the 2D hits
293  // We'll get 3D information for each of these by using the axis as a reference and use
294  // the point of closest approach as the 3D position
295 
296  // Define elements of our covariance matrix
297  float xi2(0.);
298  float xiyi(0.);
299  float xizi(0.);
300  float yi2(0.);
301  float yizi(0.);
302  float zi2(0.);
303 
304  float aveHitDoca(0.);
305  Eigen::Vector3f avePosUpdate(0., 0., 0.);
306  int nHits(0);
307 
308  // Recover existing line parameters for current cluster
309  const reco::PrincipalComponents& inputPca = pca;
310  Eigen::Vector3f avePosition(inputPca.getAvePosition());
311  Eigen::Vector3f axisDirVec(inputPca.getEigenVectors().row(0));
312 
313  // We float loop here so we can use this method for both the first time through
314  // and a second time through where we re-calculate the mean position
315  // So, we need to keep track of the poca which we do with a float vector
316  std::vector<Eigen::Vector3f> hitPosVec;
317 
318  // Outer loop over 3D hits
319  for (const auto& hit3D : hitPairVector) {
320  // Inner loop over 2D hits
321  for (const auto& hit : hit3D->getHits()) {
322  // Step one is to set up to determine the point of closest approach of this 2D hit to
323  // the cluster's current axis.
324  // Get this wire's geometry object
325  const geo::WireID& hitID = hit->WireID();
326  const geo::WireGeo& wire_geom = m_geometry->WireIDToWireGeo(hitID);
327 
328  // From this, get the parameters of the line for the wire
329  double wirePosArr[3] = {0., 0., 0.};
330  wire_geom.GetCenter(wirePosArr);
331 
332  Eigen::Vector3f wireCenter(wirePosArr[0], wirePosArr[1], wirePosArr[2]);
333  Eigen::Vector3f wireDirVec(
334  wire_geom.Direction().X(), wire_geom.Direction().Y(), wire_geom.Direction().Z());
335 
336  // Correct the wire position in x to set to correspond to the drift time
337  float hitPeak(hit->getHit()->PeakTime());
338 
339  Eigen::Vector3f wirePos(
340  detProp.ConvertTicksToX(hitPeak, hitID.Plane, hitID.TPC, hitID.Cryostat),
341  wireCenter[1],
342  wireCenter[2]);
343 
344  // Compute the wire plane normal for this view
345  Eigen::Vector3f xAxis(1., 0., 0.);
346  Eigen::Vector3f planeNormal =
347  xAxis.cross(wireDirVec); // This gives a normal vector in +z for a Y wire
348 
349  float docaInPlane(wirePos[0] - avePosition[0]);
350  float arcLenToPlane(0.);
351  float cosAxisToPlaneNormal = axisDirVec.dot(planeNormal);
352 
353  Eigen::Vector3f axisPlaneIntersection = wirePos;
354  Eigen::Vector3f hitPosTVec = wirePos;
355 
356  if (fabs(cosAxisToPlaneNormal) > 0.) {
357  Eigen::Vector3f deltaPos = wirePos - avePosition;
358 
359  arcLenToPlane = deltaPos.dot(planeNormal) / cosAxisToPlaneNormal;
360  axisPlaneIntersection = avePosition + arcLenToPlane * axisDirVec;
361  docaInPlane = wirePos[0] - axisPlaneIntersection[0];
362 
363  Eigen::Vector3f axisToInter = axisPlaneIntersection - wirePos;
364  float arcLenToDoca = axisToInter.dot(wireDirVec);
365 
366  hitPosTVec += arcLenToDoca * wireDirVec;
367  }
368 
369  // Get a vector from the wire position to our cluster's current average position
370  Eigen::Vector3f wVec = avePosition - wirePos;
371 
372  // Get the products we need to compute the arc lengths to the distance of closest approach
373  float a(axisDirVec.dot(axisDirVec));
374  float b(axisDirVec.dot(wireDirVec));
375  float c(wireDirVec.dot(wireDirVec));
376  float d(axisDirVec.dot(wVec));
377  float e(wireDirVec.dot(wVec));
378 
379  float den(a * c - b * b);
380  float arcLen1(0.);
381  float arcLen2(0.);
382 
383  // Parallel lines is a special case
384  if (fabs(den) > m_parallel) {
385  arcLen1 = (b * e - c * d) / den;
386  arcLen2 = (a * e - b * d) / den;
387  }
388  else {
389  mf::LogDebug("Cluster3D") << "** Parallel lines, need a solution here" << std::endl;
390  break;
391  }
392 
393  // Now get the hit position we'll use for the pca analysis
394  //float hitPos[] = {wirePos[0] + arcLen2 * wireDirVec[0],
395  // wirePos[1] + arcLen2 * wireDirVec[1],
396  // wirePos[2] + arcLen2 * wireDirVec[2]};
397  //float axisPos[] = {avePosition[0] + arcLen1 * axisDirVec[0],
398  // avePosition[1] + arcLen1 * axisDirVec[1],
399  // avePosition[2] + arcLen1 * axisDirVec[2]};
400  Eigen::Vector3f hitPos = wirePos + arcLen2 * wireDirVec;
401  Eigen::Vector3f axisPos = avePosition + arcLen1 * axisDirVec;
402  float deltaX = hitPos(0) - axisPos(0);
403  float deltaY = hitPos(1) - axisPos(1);
404  float deltaZ = hitPos(2) - axisPos(2);
405  float doca2 = deltaX * deltaX + deltaY * deltaY + deltaZ * deltaZ;
406  float doca = sqrt(doca2);
407 
408  docaInPlane = doca;
409 
410  aveHitDoca += fabs(docaInPlane);
411 
412  //Eigen::Vector3f deltaPos = hitPos - hitPosTVec;
413  //float deltaDoca = doca - docaInPlane;
414 
415  //if (fabs(deltaPos[0]) > 1. || fabs(deltaPos[1]) > 1. || fabs(deltaPos[2]) > 1. || fabs(deltaDoca) > 2.)
416  //{
417  // std::cout << "**************************************************************************************" << std::endl;
418  // std::cout << "Diff in doca: " << deltaPos[0] << "," << deltaPos[1] << "," << deltaPos[2] << ", deltaDoca: " << deltaDoca << std::endl;
419  // std::cout << "-- HitPosTVec: " << hitPosTVec[0] << "," << hitPosTVec[1] << "," << hitPosTVec[2] << std::endl;
420  // std::cout << "-- hitPos: " << hitPos[0] << "," << hitPos[1] << "," << hitPos[2] << std::endl;
421  // std::cout << "-- WirePos: " << wirePos[0] << "," << wirePos[1] << "," << wirePos[2] << std::endl;
422  //}
423 
424  // Set the hit's doca and arclen
425  hit->setDocaToAxis(fabs(docaInPlane));
426  hit->setArcLenToPoca(arcLenToPlane);
427 
428  // If this point is considered an outlier then we skip
429  // the accumulation for average position and covariance
430  if (hit->getStatusBits() & 0x80) { continue; }
431 
432  //hitPosTVec = hitPos;
433 
434  avePosUpdate += hitPosTVec;
435 
436  hitPosVec.push_back(hitPosTVec);
437 
438  nHits++;
439  }
440  }
441 
442  // Get updated average position
443  avePosUpdate /= float(nHits);
444 
445  // Get the average hit doca
446  aveHitDoca /= float(nHits);
447 
448  if (updateAvePos) { avePosition = avePosUpdate; }
449 
450  // Now loop through the hits and build out the covariance matrix
451  for (auto& hitPos : hitPosVec) {
452  // And increment the values in the covariance matrix
453  float x = hitPos[0] - avePosition[0];
454  float y = hitPos[1] - avePosition[1];
455  float z = hitPos[2] - avePosition[2];
456 
457  xi2 += x * x;
458  xiyi += x * y;
459  xizi += x * z;
460  yi2 += y * y;
461  yizi += y * z;
462  zi2 += z * z;
463  }
464 
465  // Accumulation done, now do the actual work
466 
467  // Using Eigen package
468  Eigen::Matrix3f sig;
469 
470  sig << xi2, xiyi, xizi, xiyi, yi2, yizi, xizi, yizi, zi2;
471 
472  sig *= 1. / (nHits - 1);
473 
474  Eigen::SelfAdjointEigenSolver<Eigen::Matrix3f> eigenMat(sig);
475 
476  if (eigenMat.info() == Eigen::ComputationInfo::Success) {
477  // Now copy output
478  // The returned eigen values and vectors will be returned in an xyz system where x is the smallest spread,
479  // y is the next smallest and z is the largest. Adopt that convention going forward
480  reco::PrincipalComponents::EigenValues recobEigenVals = eigenMat.eigenvalues().cast<float>();
482  eigenMat.eigenvectors().transpose().cast<float>();
483 
484  // Store away
486  true, nHits, recobEigenVals, recobEigenVecs, avePosition, aveHitDoca);
487  }
488  else {
489  mf::LogDebug("Cluster3D") << "PCA decompose failure, numPairs = " << nHits << std::endl;
491  }
492 
493  return;
494  }
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:65
process_name opflash particleana ie ie ie z
process_name opflash particleana ie x
CryostatID_t Cryostat
Index of cryostat.
Definition: geo_types.h:212
process_name hit
Definition: cheaterreco.fcl:51
process_name gaushit a
process_name opflash particleana ie ie y
Eigen::Vector3f EigenValues
Definition: Cluster3D.h:225
const Eigen::Vector3f & getAvePosition() const
Definition: Cluster3D.h:247
PlaneID_t Plane
Index of the plane within its TPC.
Definition: geo_types.h:493
Vector Direction() const
Definition: WireGeo.h:587
double ConvertTicksToX(double ticks, int p, int t, int c) const
Eigen::Matrix3f EigenVectors
Definition: Cluster3D.h:226
do i e
void GetCenter(double *xyz, double localz=0.0) const
Fills the world coordinate of a point on the wire.
Definition: WireGeo.cxx:73
TPCID_t TPC
Index of the TPC within its cryostat.
Definition: geo_types.h:406
const EigenVectors & getEigenVectors() const
Definition: Cluster3D.h:246
WireGeo const & WireIDToWireGeo(geo::WireID const &wireid) const
void lar_cluster3d::PrincipalComponentsAlg::PCAAnalysis_3D ( const reco::HitPairListPtr hitPairList,
reco::PrincipalComponents pca,
bool  skeletonOnly = false 
) const

Definition at line 141 of file PrincipalComponentsAlg.cxx.

144  {
145  // We want to run a PCA on the input TkrVecPoints...
146  // The steps are:
147  // 1) do a mean normalization of the input vec points
148  // 2) compute the covariance matrix
149  // 3) run the SVD
150  // 4) extract the eigen vectors and values
151  // see what happens
152 
153  // Run through the HitPairList and get the mean position of all the hits
154  Eigen::Vector3d meanPos(Eigen::Vector3d::Zero());
155  double meanWeightSum(0.);
156  int numPairsInt(0);
157 
158  // const float minimumDeltaPeakSig(0.00001);
159  double minimumDeltaPeakSig(0.00001);
160 
161  // Want to use the hit "chi square" to weight the hits but we need to put a lower limit on its value
162  // to prevent a few hits being over counted.
163  // This is a bit experimental until we can evaluate the cost (time to calculate) vs the benefit
164  // (better fits)..
165  std::vector<double> hitChiSquareVec;
166 
167  hitChiSquareVec.resize(hitPairVector.size());
168 
169  std::transform(hitPairVector.begin(),
170  hitPairVector.end(),
171  hitChiSquareVec.begin(),
172  [](const auto& hit) { return hit->getHitChiSquare(); });
173  std::sort(hitChiSquareVec.begin(), hitChiSquareVec.end());
174 
175  size_t numToKeep = 0.8 * hitChiSquareVec.size();
176 
177  hitChiSquareVec.resize(numToKeep);
178 
179  double aveValue = std::accumulate(hitChiSquareVec.begin(), hitChiSquareVec.end(), double(0.)) /
180  double(hitChiSquareVec.size());
181  double rms = std::sqrt(std::inner_product(hitChiSquareVec.begin(),
182  hitChiSquareVec.end(),
183  hitChiSquareVec.begin(),
184  0.,
185  std::plus<>(),
186  [aveValue](const auto& left, const auto& right) {
187  return (left - aveValue) * (right - aveValue);
188  }) /
189  double(hitChiSquareVec.size()));
190 
191  minimumDeltaPeakSig = std::max(minimumDeltaPeakSig, aveValue - rms);
192 
193  // std::cout << "===>> Calculating PCA, ave chiSquare: " << aveValue << ", rms: " << rms << ", cut: " << minimumDeltaPeakSig << std::endl;
194 
195  for (const auto& hit : hitPairVector) {
196  if (skeletonOnly && !((hit->getStatusBits() & reco::ClusterHit3D::SKELETONHIT) ==
198  continue;
199 
200  // Weight the hit by the peak time difference significance
201  double weight = std::max(minimumDeltaPeakSig,
202  double(hit->getHitChiSquare())); // hit->getDeltaPeakTime());
203  // ///hit->getSigmaPeakTime());
204 
205  meanPos(0) += hit->getPosition()[0] * weight;
206  meanPos(1) += hit->getPosition()[1] * weight;
207  meanPos(2) += hit->getPosition()[2] * weight;
208  numPairsInt++;
209 
210  meanWeightSum += weight;
211  }
212 
213  meanPos /= meanWeightSum;
214 
215  // Define elements of our covariance matrix
216  double xi2(0.);
217  double xiyi(0.);
218  double xizi(0.0);
219  double yi2(0.0);
220  double yizi(0.0);
221  double zi2(0.);
222  double weightSum(0.);
223 
224  // Back through the hits to build the matrix
225  for (const auto& hit : hitPairVector) {
226  if (skeletonOnly && !((hit->getStatusBits() & reco::ClusterHit3D::SKELETONHIT) ==
228  continue;
229 
230  double weight = 1. / std::max(minimumDeltaPeakSig,
231  double(hit->getHitChiSquare())); // hit->getDeltaPeakTime());
232  // ///hit->getSigmaPeakTime());
233  double x = (hit->getPosition()[0] - meanPos(0)) * weight;
234  double y = (hit->getPosition()[1] - meanPos(1)) * weight;
235  double z = (hit->getPosition()[2] - meanPos(2)) * weight;
236 
237  weightSum += weight * weight;
238 
239  xi2 += x * x;
240  xiyi += x * y;
241  xizi += x * z;
242  yi2 += y * y;
243  yizi += y * z;
244  zi2 += z * z;
245  }
246 
247  // Using Eigen package
248  Eigen::Matrix3d sig;
249 
250  sig << xi2, xiyi, xizi, xiyi, yi2, yizi, xizi, yizi, zi2;
251 
252  sig *= 1. / weightSum;
253 
254  Eigen::SelfAdjointEigenSolver<Eigen::Matrix3d> eigenMat(sig);
255 
256  if (eigenMat.info() == Eigen::ComputationInfo::Success) {
257  // Now copy output
258  // The returned eigen values and vectors will be returned in an xyz system where x is the smallest spread,
259  // y is the next smallest and z is the largest. Adopt that convention going forward
260  reco::PrincipalComponents::EigenValues recobEigenVals = eigenMat.eigenvalues().cast<float>();
262  eigenMat.eigenvectors().transpose().cast<float>();
263 
264  // Check for a special case (which may have gone away with switch back to doubles for computation?)
265  if (std::isnan(recobEigenVals[0])) {
266  std::cout << "==> Third eigenvalue returns a nan" << std::endl;
267 
268  recobEigenVals[0] = 0.;
269 
270  // Assume the third axis is also kaput?
271  recobEigenVecs.row(0) = recobEigenVecs.row(1).cross(recobEigenVecs.row(2));
272  }
273 
274  // Store away
276  true, numPairsInt, recobEigenVals, recobEigenVecs, meanPos.cast<float>());
277  }
278  else {
279  mf::LogDebug("Cluster3D") << "PCA decompose failure, numPairs = " << numPairsInt << std::endl;
281  }
282 
283  return;
284  }
process_name opflash particleana ie ie ie z
static constexpr Sample_t transform(Sample_t sample)
process_name opflash particleana ie x
walls no right
Definition: selectors.fcl:105
process_name hit
Definition: cheaterreco.fcl:51
process_name opflash particleana ie ie y
Eigen::Vector3f EigenValues
Definition: Cluster3D.h:225
walls no left
Definition: selectors.fcl:105
Eigen::Matrix3f EigenVectors
Definition: Cluster3D.h:226
Hit is a &quot;skeleton&quot; hit.
Definition: Cluster3D.h:99
BEGIN_PROLOG could also be cout
void lar_cluster3d::PrincipalComponentsAlg::PCAAnalysis_calc2DDocas ( const reco::Hit2DListPtr hit2DVector,
const reco::PrincipalComponents pca 
) const

Definition at line 552 of file PrincipalComponentsAlg.cxx.

554  {
555  // Our mission, should we choose to accept it, is to scan through the 2D hits and reject
556  // any outliers. Basically, any hit outside a scaled range of the average doca from the
557  // first pass is marked by setting the bit in the status word.
558 
559  // We'll need the current PCA axis to determine doca and arclen
560  Eigen::Vector3f avePosition(
561  pca.getAvePosition()[0], pca.getAvePosition()[1], pca.getAvePosition()[2]);
562  Eigen::Vector3f axisDirVec(pca.getEigenVectors().row(2));
563 
564  // Recover the principle eigen value for range constraints
565  float maxArcLen = 4. * sqrt(pca.getEigenValues()[0]);
566 
567  // We want to keep track of the average
568  float aveHitDoca(0.);
569 
570  // Outer loop over views
571  for (const auto* hit : hit2DListPtr) {
572  // Step one is to set up to determine the point of closest approach of this 2D hit to
573  // the cluster's current axis. We do that by finding the point of intersection of the
574  // cluster's axis with a plane defined by the wire the hit is associated with.
575  // Get this wire's geometry object
576  const geo::WireID& hitID = hit->WireID();
577  const geo::WireGeo& wire_geom = m_geometry->WireIDToWireGeo(hitID);
578 
579  // From this, get the parameters of the line for the wire
580  double wirePosArr[3] = {0., 0., 0.};
581  wire_geom.GetCenter(wirePosArr);
582 
583  Eigen::Vector3f wireCenter(wirePosArr[0], wirePosArr[1], wirePosArr[2]);
584  Eigen::Vector3f wireDirVec(
585  wire_geom.Direction().X(), wire_geom.Direction().Y(), wire_geom.Direction().Z());
586 
587  // Correct the wire position in x to set to correspond to the drift time
588  Eigen::Vector3f wirePos(hit->getXPosition(), wireCenter[1], wireCenter[2]);
589 
590  // Compute the wire plane normal for this view
591  Eigen::Vector3f xAxis(1., 0., 0.);
592  Eigen::Vector3f planeNormal =
593  xAxis.cross(wireDirVec); // This gives a normal vector in +z for a Y wire
594 
595  float arcLenToPlane(0.);
596  float docaInPlane(wirePos[0] - avePosition[0]);
597  float cosAxisToPlaneNormal = axisDirVec.dot(planeNormal);
598 
599  Eigen::Vector3f axisPlaneIntersection = wirePos;
600 
601  // If current cluster axis is not parallel to wire plane then find intersection point
602  if (fabs(cosAxisToPlaneNormal) > 0.) {
603  Eigen::Vector3f deltaPos = wirePos - avePosition;
604 
605  arcLenToPlane =
606  std::min(float(deltaPos.dot(planeNormal) / cosAxisToPlaneNormal), maxArcLen);
607  axisPlaneIntersection = avePosition + arcLenToPlane * axisDirVec;
608 
609  Eigen::Vector3f axisToInter = axisPlaneIntersection - wirePos;
610  float arcLenToDoca = axisToInter.dot(wireDirVec);
611 
612  // If the arc length along the wire to the poca is outside the TPC then reset
613  if (fabs(arcLenToDoca) > wire_geom.HalfL()) arcLenToDoca = wire_geom.HalfL();
614 
615  // If we were successful in getting to the wire plane then the doca is simply the
616  // difference in x coordinates... but we hvae to worry about the special cases so
617  // we calculate a 3D doca based on arclengths above...
618  Eigen::Vector3f docaVec = axisPlaneIntersection - (wirePos + arcLenToDoca * wireDirVec);
619  docaInPlane = docaVec.norm();
620  }
621 
622  aveHitDoca += fabs(docaInPlane);
623 
624  // Set the hit's doca and arclen
625  hit->setDocaToAxis(fabs(docaInPlane));
626  hit->setArcLenToPoca(arcLenToPlane);
627  }
628 
629  // Compute the average and store
630  aveHitDoca /= float(hit2DListPtr.size());
631 
632  pca.setAveHitDoca(aveHitDoca);
633 
634  return;
635  }
Geometry description of a TPC wireThe wire is a single straight segment on a wire plane...
Definition: WireGeo.h:65
void setAveHitDoca(double doca) const
Definition: Cluster3D.h:251
process_name hit
Definition: cheaterreco.fcl:51
const EigenValues & getEigenValues() const
Definition: Cluster3D.h:245
const Eigen::Vector3f & getAvePosition() const
Definition: Cluster3D.h:247
Vector Direction() const
Definition: WireGeo.h:587
double HalfL() const
Returns half the length of the wire [cm].
Definition: WireGeo.h:128
void GetCenter(double *xyz, double localz=0.0) const
Fills the world coordinate of a point on the wire.
Definition: WireGeo.cxx:73
const EigenVectors & getEigenVectors() const
Definition: Cluster3D.h:246
WireGeo const & WireIDToWireGeo(geo::WireID const &wireid) const
void lar_cluster3d::PrincipalComponentsAlg::PCAAnalysis_calc3DDocas ( const reco::HitPairListPtr hitPairVector,
const reco::PrincipalComponents pca 
) const

Definition at line 497 of file PrincipalComponentsAlg.cxx.

499  {
500  // Our mission, should we choose to accept it, is to scan through the 2D hits and reject
501  // any outliers. Basically, any hit outside a scaled range of the average doca from the
502  // first pass is marked by setting the bit in the status word.
503 
504  // We'll need the current PCA axis to determine doca and arclen
505  Eigen::Vector3f avePosition(
506  pca.getAvePosition()[0], pca.getAvePosition()[1], pca.getAvePosition()[2]);
507  Eigen::Vector3f axisDirVec(pca.getEigenVectors().row(2));
508 
509  // We want to keep track of the average
510  float aveDoca3D(0.);
511 
512  // Outer loop over views
513  for (const auto* clusterHit3D : hitPairVector) {
514  // Always reset the existing status bit
515  clusterHit3D->clearStatusBits(0x80);
516 
517  // Now we need to calculate the doca and poca...
518  // Start by getting this hits position
519  Eigen::Vector3f clusPos(clusterHit3D->getPosition()[0],
520  clusterHit3D->getPosition()[1],
521  clusterHit3D->getPosition()[2]);
522 
523  // Form a TVector from this to the cluster average position
524  Eigen::Vector3f clusToHitVec = clusPos - avePosition;
525 
526  // With this we can get the arclength to the doca point
527  float arclenToPoca = clusToHitVec.dot(axisDirVec);
528 
529  // Get the coordinates along the axis for this point
530  Eigen::Vector3f docaPos = avePosition + arclenToPoca * axisDirVec;
531 
532  // Now get doca and poca
533  Eigen::Vector3f docaPosToClusPos = clusPos - docaPos;
534  float docaToAxis = docaPosToClusPos.norm();
535 
536  aveDoca3D += docaToAxis;
537 
538  // Ok, set the values in the hit
539  clusterHit3D->setDocaToAxis(docaToAxis);
540  clusterHit3D->setArclenToPoca(arclenToPoca);
541  }
542 
543  // Compute the average and store
544  aveDoca3D /= float(hitPairVector.size());
545 
546  pca.setAveHitDoca(aveDoca3D);
547 
548  return;
549  }
void setAveHitDoca(double doca) const
Definition: Cluster3D.h:251
const Eigen::Vector3f & getAvePosition() const
Definition: Cluster3D.h:247
const EigenVectors & getEigenVectors() const
Definition: Cluster3D.h:246
int lar_cluster3d::PrincipalComponentsAlg::PCAAnalysis_reject2DOutliers ( const reco::HitPairListPtr hitPairVector,
reco::PrincipalComponents pca,
float  aveHitDoca 
) const

Definition at line 638 of file PrincipalComponentsAlg.cxx.

641  {
642  // Our mission, should we choose to accept it, is to scan through the 2D hits and reject
643  // any outliers. Basically, any hit outside a scaled range of the average doca from the
644  // first pass is marked by setting the bit in the status word.
645 
646  // First get the average doca scaled by some appropriate factor
647  int numRejHits(0);
648 
649  // Outer loop over views
650  for (const auto& hit3D : hitPairVector) {
651  // Inner loop over hits in this view
652  for (const auto& hit : hit3D->getHits()) {
653  // Always reset the existing status bit
654  hit->clearStatusBits(0x80);
655 
656  if (hit->getDocaToAxis() > maxDocaAllowed) {
657  hit->setStatusBit(0x80);
658  numRejHits++;
659  }
660  }
661  }
662 
663  return numRejHits;
664  }
process_name hit
Definition: cheaterreco.fcl:51
int lar_cluster3d::PrincipalComponentsAlg::PCAAnalysis_reject3DOutliers ( const reco::HitPairListPtr hitPairVector,
const reco::PrincipalComponents pca,
float  aveHitDoca 
) const

Definition at line 667 of file PrincipalComponentsAlg.cxx.

670  {
671  // Our mission, should we choose to accept it, is to scan through the 2D hits and reject
672  // any outliers. Basically, any hit outside a scaled range of the average doca from the
673  // first pass is marked by setting the bit in the status word.
674 
675  // First get the average doca scaled by some appropriate factor
676  int numRejHits(0);
677 
678  // We'll need the current PCA axis to determine doca and arclen
679  Eigen::Vector3f avePosition(
680  pca.getAvePosition()[0], pca.getAvePosition()[1], pca.getAvePosition()[2]);
681  Eigen::Vector3f axisDirVec(pca.getEigenVectors().row(2));
682 
683  // Outer loop over views
684  for (const auto* clusterHit3D : hitPairVector) {
685  // Always reset the existing status bit
686  clusterHit3D->clearStatusBits(0x80);
687 
688  // Now we need to calculate the doca and poca...
689  // Start by getting this hits position
690  Eigen::Vector3f clusPos(clusterHit3D->getPosition()[0],
691  clusterHit3D->getPosition()[1],
692  clusterHit3D->getPosition()[2]);
693 
694  // Form a TVector from this to the cluster average position
695  Eigen::Vector3f clusToHitVec = clusPos - avePosition;
696 
697  // With this we can get the arclength to the doca point
698  float arclenToPoca = clusToHitVec.dot(axisDirVec);
699 
700  // Get the coordinates along the axis for this point
701  Eigen::Vector3f docaPos = avePosition + arclenToPoca * axisDirVec;
702 
703  // Now get doca and poca
704  Eigen::Vector3f docaPosToClusPos = clusPos - docaPos;
705  float docaToAxis = docaPosToClusPos.norm();
706 
707  // Ok, set the values in the hit
708  clusterHit3D->setDocaToAxis(docaToAxis);
709  clusterHit3D->setArclenToPoca(arclenToPoca);
710 
711  // Check to see if this is a keeper
712  if (clusterHit3D->getDocaToAxis() > maxDocaAllowed) {
713  clusterHit3D->setStatusBit(0x80);
714  numRejHits++;
715  }
716  }
717 
718  return numRejHits;
719  }
const Eigen::Vector3f & getAvePosition() const
Definition: Cluster3D.h:247
const EigenVectors & getEigenVectors() const
Definition: Cluster3D.h:246

Member Data Documentation

const geo::Geometry* lar_cluster3d::PrincipalComponentsAlg::m_geometry
private

Definition at line 75 of file PrincipalComponentsAlg.h.

float lar_cluster3d::PrincipalComponentsAlg::m_parallel
private

means lines are parallel

Definition at line 74 of file PrincipalComponentsAlg.h.


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