All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Cluster3D.h
Go to the documentation of this file.
1 /**
2  *
3  * @brief Definition of utility objects for use in the 3D clustering for LArSoft
4  *
5  * The objects defined in this module are intended for internal use by the
6  * 3D clustering (see Cluster3D_module.cc in larreco/ClusterFinder).
7  * These objects mostly contain volatile information and are not suitable
8  * for storage in the art event store
9  *
10  * @author usher@slac.stanford.edu
11  *
12  */
13 
14 #ifndef RECO_CLUSTER3D_H
15 #define RECO_CLUSTER3D_H
16 
17 #include <iosfwd>
18 #include <vector>
19 #include <list>
20 #include <set>
21 #include <map>
22 #include <unordered_map>
23 
26 namespace recob { class Hit; }
27 
28 // Eigen
29 #include <Eigen/Core>
30 
31 namespace reco {
32 
33 // First define a container object to augment the sparse 2D hit information
35 {
36 public:
37 
38  ClusterHit2D(); // Default constructor
39 
40 private:
41 
42  mutable unsigned m_statusBits; ///< Volatile status information of this 3D hit
43  mutable float m_docaToAxis; ///< DOCA of hit at POCA to associated cluster axis
44  mutable float m_arcLenToPoca; ///< arc length to POCA along cluster axis
45  float m_xPosition; ///< The x coordinate for this hit
46  float m_timeTicks; ///< The time (in ticks) for this hit
47  geo::WireID m_wireID; ///< Keep track this particular hit's wireID
48  const recob::Hit* m_hit; ///< Hit we are augmenting
49 
50 public:
51 
52  enum StatusBits { SHAREDINPAIR = 0x00080000,
53  SHAREDINTRIPLET = 0x00040000,
54  USEDINPAIR = 0x00008000,
55  USEDINTRIPLET = 0x00004000,
56  SHAREDINCLUSTER = 0x00000200,
57  USEDINCLUSTER = 0x00000100,
58  USED = 0x00000001
59  };
60 
61  ClusterHit2D(unsigned statusBits,
62  float doca,
63  float poca,
64  float xPosition,
65  float timeTicks,
66  const geo::WireID& wireID,
67  const recob::Hit* recobHit);
68 
69  ClusterHit2D(const ClusterHit2D&);
70 
71  unsigned getStatusBits() const {return m_statusBits;}
72  float getDocaToAxis() const {return m_docaToAxis;}
73  float getArcLenToPoca() const {return m_arcLenToPoca;}
74  float getXPosition() const {return m_xPosition;}
75  float getTimeTicks() const {return m_timeTicks;}
76  const geo::WireID& WireID() const {return m_wireID;}
77  const recob::Hit* getHit() const {return m_hit;}
78 
79  void setStatusBit(unsigned bits) const {m_statusBits |= bits;}
80  void clearStatusBits(unsigned bits) const {m_statusBits &= ~bits;}
81  void setDocaToAxis(float doca) const {m_docaToAxis = doca;}
82  void setArcLenToPoca(float poca) const {m_arcLenToPoca = poca;}
83 
84  void setHit(const recob::Hit* hit) {m_hit = hit;}
85 
86  friend std::ostream& operator << (std::ostream& o, const ClusterHit2D& c);
87  friend bool operator < (const ClusterHit2D & a, const ClusterHit2D & b);
88 
89 };
90 
91 using ClusterHit2DVec = std::vector<const reco::ClusterHit2D*>;
92 
93 // Now define an object with the recob::Hit information that will comprise the 3D cluster
95 {
96 public:
97 
98  enum StatusBits { REJECTEDHIT = 0x80000000, ///< Hit has been rejected for any reason
99  SKELETONHIT = 0x10000000, ///< Hit is a "skeleton" hit
100  EDGEHIT = 0x20000000, ///< Hit is an "edge" hit
101  SEEDHIT = 0x40000000, ///< Hit is part of Seed for track fits
102  MADESPACEPOINT = 0x08000000, ///< Hit has been made into Space Point
103  CONVEXHULLVTX = 0x04000000, ///< Point is on primary cluster convex hull
104  EXTREMEPOINT = 0x02000000, ///< Is a convex hull extreme point
105  SKELETONPOSAVE = 0x00100000, ///< Skeleton hit position averaged
106  CLUSTERVISITED = 0x00008000, ///< "visited" by a clustering algorithm
107  CLUSTERNOISE = 0x00004000, ///< Labelled "noise" by a clustering algorithm
108  CLUSTERATTACHED = 0x00002000, ///< attached to a cluster
109  CLUSTERSHARED = 0x00001000, ///< 3D hit has 2D hits shared between clusters
110  PATHCHECKED = 0x00000800, ///< Path checking algorithm has seen this hit
111  SELECTEDBYMST = 0x00000100, ///< Hit has been used in Cluster Splitting MST
112  PCAOUTLIER = 0x00000080, ///< Hit labelled outlier in PCA
113  HITINVIEW0 = 0x00000001, ///< Hit contains 2D hit from view 0 (u plane)
114  HITINVIEW1 = 0x00000002, ///< Hit contains 2D hit from view 1 (v plane)
115  HITINVIEW2 = 0x00000004 ///< Hit contains 2D hit from view 2 (w plane)
116  };
117 
118 
119  ClusterHit3D(); // Default constructor
120 
121  ClusterHit3D(size_t id,
122  unsigned int statusBits,
123  const Eigen::Vector3f& position,
124  float totalCharge,
125  float avePeakTime,
126  float deltaPeakTime,
127  float sigmaPeakTime,
128  float hitChiSquare,
129  float overlapFraction,
130  float chargeAsymmetry,
131  float docaToAxis,
132  float arclenToPoca,
133  const ClusterHit2DVec& hitVec,
134  const std::vector<float>& hitDelTSigVec,
135  const std::vector<geo::WireID>& wireIDVec);
136 
137  ClusterHit3D(const ClusterHit3D&);
138 
139  void initialize(size_t id,
140  unsigned int statusBits,
141  const Eigen::Vector3f& position,
142  float totalCharge,
143  float avePeakTime,
144  float deltaPeakTime,
145  float sigmaPeakTime,
146  float hitChiSquare,
147  float overlapFraction,
148  float chargeAsymmetry,
149  float docaToAxis,
150  float arclenToPoca,
151  const ClusterHit2DVec& hitVec,
152  const std::vector<float>& hitDelTSigVec,
153  const std::vector<geo::WireID>& wireIDVec);
154 
155  size_t getID() const {return fID;}
156  unsigned int getStatusBits() const {return fStatusBits;}
157  const Eigen::Vector3f getPosition() const {return fPosition;}
158  float getX() const {return fPosition[0];}
159  float getY() const {return fPosition[1];}
160  float getZ() const {return fPosition[2];}
161  float getTotalCharge() const {return fTotalCharge;}
162  float getAvePeakTime() const {return fAvePeakTime;}
163  float getDeltaPeakTime() const {return fDeltaPeakTime;}
164  float getSigmaPeakTime() const {return fSigmaPeakTime;}
165  float getHitChiSquare() const {return fHitChiSquare;}
166  float getOverlapFraction() const {return fOverlapFraction;}
167  float getChargeAsymmetry() const {return fChargeAsymmetry;}
168  float getDocaToAxis() const {return fDocaToAxis;}
169  float getArclenToPoca() const {return fArclenToPoca;}
170  const ClusterHit2DVec& getHits() const {return fHitVector;}
171  const std::vector<float> getHitDelTSigVec() const {return fHitDelTSigVec;}
172  const std::vector<geo::WireID>& getWireIDs() const {return fWireIDVector;}
173 
175 
176  bool bitsAreSet(const unsigned int& bitsToCheck) const {return fStatusBits & bitsToCheck;}
177 
178  void setID(const size_t& id) const {fID = id;}
179  void setStatusBit(unsigned bits) const {fStatusBits |= bits;}
180  void clearStatusBits(unsigned bits) const {fStatusBits &= ~bits;}
181  void setDocaToAxis(double doca) const {fDocaToAxis = doca;}
182  void setArclenToPoca(double poca) const {fArclenToPoca = poca;}
183  void setWireID(const geo::WireID& wid) const;
184 
185  void setPosition(const Eigen::Vector3f& pos) const {fPosition = pos;}
186 
187  const bool operator<(const reco::ClusterHit3D& other) const
188  {
189  if (fPosition[2] != other.fPosition[2]) return fPosition[2] < other.fPosition[2];
190  else return fPosition[0] < other.fPosition[0];
191  }
192 
193  const bool operator==(const reco::ClusterHit3D& other) const
194  {
195  return fID == other.fID;
196  }
197 
198  friend std::ostream& operator << (std::ostream& o, const ClusterHit3D& c);
199  //friend bool operator < (const ClusterHit3D & a, const ClusterHit3D & b);
200 
201 private:
202 
203  mutable size_t fID; ///< "id" of this hit (useful for indexing)
204  mutable unsigned int fStatusBits; ///< Volatile status information of this 3D hit
205  mutable Eigen::Vector3f fPosition; ///< position of this hit combination in world coordinates
206  float fTotalCharge; ///< Sum of charges of all associated recob::Hits
207  float fAvePeakTime; ///< Average peak time of all associated recob::Hits
208  float fDeltaPeakTime; ///< Largest delta peak time of associated recob::Hits
209  float fSigmaPeakTime; ///< Quad sum of peak time sigmas
210  float fHitChiSquare; ///< Hit ChiSquare relative to the average time
211  float fOverlapFraction; ///< Hit overlap fraction start/stop of triplet
212  float fChargeAsymmetry; ///< Assymetry of average of two closest to third charge
213  mutable float fDocaToAxis; ///< DOCA to the associated cluster axis
214  mutable float fArclenToPoca; ///< arc length along axis to DOCA point
215  ClusterHit2DVec fHitVector; ///< Hits comprising this 3D hit
216  mutable std::vector<float> fHitDelTSigVec; ///< Delta t of hit to matching pair / sig
217  mutable std::vector<geo::WireID> fWireIDVector; ///< Wire ID's for the planes making up hit
218 };
219 
220 // We also need to define a container for the output of the PCA Analysis
222 {
223 public:
224 
225  using EigenValues = Eigen::Vector3f;
226  using EigenVectors = Eigen::Matrix3f;
227 
229 
230 private:
231 
232  bool m_svdOK; ///< SVD Decomposition was successful
233  int m_numHitsUsed; ///< Number of hits in the decomposition
234  EigenValues m_eigenValues; ///< Eigen values from SVD decomposition
235  EigenVectors m_eigenVectors; ///< The three principle axes
236  Eigen::Vector3f m_avePosition; ///< Average position of hits fed to PCA
237  mutable double m_aveHitDoca; ///< Average doca of hits used in PCA
238 
239 public:
240 
241  PrincipalComponents(bool ok, int nHits, const EigenValues& eigenValues, const EigenVectors& eigenVecs, const Eigen::Vector3f& avePos, const float aveHitDoca = 9999.);
242 
243  bool getSvdOK() const {return m_svdOK;}
244  int getNumHitsUsed() const {return m_numHitsUsed;}
245  const EigenValues& getEigenValues() const {return m_eigenValues;}
246  const EigenVectors& getEigenVectors() const {return m_eigenVectors;}
247  const Eigen::Vector3f& getAvePosition() const {return m_avePosition;}
248  const float getAveHitDoca() const {return m_aveHitDoca;}
249 
250  void flipAxis(size_t axis);
251  void setAveHitDoca(double doca) const {m_aveHitDoca = doca;}
252 
253  friend std::ostream& operator << (std::ostream & o, const PrincipalComponents& a);
254  friend bool operator < (const PrincipalComponents& a, const PrincipalComponents& b);
255 
256 };
257 
259 {
260 public:
261 
262  Cluster3D(); ///Default constructor
263 
264 private:
265 
266  mutable unsigned m_statusBits; ///< Status bits for the cluster
267  PrincipalComponents m_pcaResults; ///< Output of the prinicipal componenets analysis
268  float m_totalCharge; ///< Total charge in the cluster
269  float m_startPosition[3]; ///< "start" position for cluster (world coordinates)
270  float m_endPosition[3]; ///< "end" position for cluster
271  int m_clusterIdx; ///< ID for this cluster
272 
273 public:
274  Cluster3D(unsigned statusBits,
275  const PrincipalComponents& pcaResults,
276  float totalCharge,
277  const float* startPosition,
278  const float* endPosition,
279  int idx);
280 
281  unsigned getStatusBits() const {return m_statusBits;}
283  float getTotalCharge() const {return m_totalCharge;}
284  const float* getStartPosition() const {return m_startPosition;}
285  const float* getEndPosition() const {return m_endPosition;}
286  int getClusterIdx() const {return m_clusterIdx;}
287 
288  void setStatusBit(unsigned bits) const {m_statusBits |= bits;}
289  void clearStatusBits(unsigned bits) const {m_statusBits &= ~bits;}
290 
292  friend std::ostream& operator << (std::ostream& o, const Cluster3D& c);
293  friend bool operator < (const Cluster3D & a, const Cluster3D & b);
294 };
295 
296 /**
297  * @brief A utility class used in construction of 3D clusters
298  */
300 {
301 public:
303  m_sigmaStartTime(1.),
304  m_endTime(0.),
305  m_sigmaEndTime(1.),
306  m_totalCharge(0.),
307  m_startWire(9999999),
308  m_endWire(0),
309  m_plane(geo::PlaneID()),
310  m_view(geo::kUnknown)
311  {
312  m_hitVector.clear();
313  }
314 
316 
317  float m_startTime;
319  float m_endTime;
322  unsigned int m_startWire;
323  unsigned int m_endWire;
327 };
328 
329 
330 /**
331  * @brief export some data structure definitions
332  */
333 using Hit2DListPtr = std::list<const reco::ClusterHit2D*>;
334 using HitPairListPtr = std::list<const reco::ClusterHit3D*>;
335 using HitPairSetPtr = std::set<const reco::ClusterHit3D*>;
336 using HitPairListPtrList = std::list<HitPairListPtr>;
337 using HitPairClusterMap = std::map<int, HitPairListPtr>;
338 using HitPairList = std::list<reco::ClusterHit3D>;
339 //using HitPairList = std::list<std::unique_ptr<reco::ClusterHit3D>>;
340 
341 using PCAHitPairClusterMapPair = std::pair<reco::PrincipalComponents, reco::HitPairClusterMap::iterator>;
342 using PlaneToClusterParamsMap = std::map<size_t, RecobClusterParameters>;
343 using EdgeTuple = std::tuple<const reco::ClusterHit3D*,const reco::ClusterHit3D*,double>;
344 using EdgeList = std::list<EdgeTuple>;
345 using Hit3DToEdgePair = std::pair<const reco::ClusterHit3D*, reco::EdgeList>;
346 using Hit3DToEdgeMap = std::unordered_map<const reco::ClusterHit3D*, reco::EdgeList>;
347 using Hit2DToHit3DListMap = std::unordered_map<const reco::ClusterHit2D*, reco::HitPairListPtr>;
348 //using VertexPoint = Eigen::Vector3f;
349 //using VertexPointList = std::list<Eigen::Vector3f>;
350 
351 using ProjectedPoint = std::tuple<float, float, const reco::ClusterHit3D*>; ///< Projected coordinates and pointer to hit
352 using ProjectedPointList = std::list<ProjectedPoint>;
353 using ConvexHullKinkTuple = std::tuple<ProjectedPoint, Eigen::Vector2f, Eigen::Vector2f>; ///< Point plus edges that point to it
354 using ConvexHullKinkTupleList = std::list<ConvexHullKinkTuple>;
355 
356 /**
357  * @brief Define a container for working with the convex hull
358  */
360 {
361 public:
363  {
364  fProjectedPointList.clear(),
365  fConvexHullPointList.clear(),
366  fConvexHullEdgeMap.clear(),
367  fConvexHullEdgeList.clear(),
368  fConvexHullExtremePoints.clear(),
369  fConvexHullKinkPoints.clear();
370  }
371 
372  void clear()
373  {
374  fProjectedPointList.clear(),
375  fConvexHullPointList.clear(),
376  fConvexHullEdgeMap.clear(),
377  fConvexHullEdgeList.clear(),
378  fConvexHullExtremePoints.clear(),
379  fConvexHullKinkPoints.clear();
380  }
381 
388 
389 
390 private:
391  reco::ProjectedPointList fProjectedPointList; ///< The input set of points projected onto plane encompassed by the hull
392  reco::ProjectedPointList fConvexHullPointList; ///< The points on the convex hull
393  reco::Hit3DToEdgeMap fConvexHullEdgeMap; ///< Map from 3D hit to associated edge
394  reco::EdgeList fConvexHullEdgeList; ///< An edge list translated back to 3D hits
395  reco::ProjectedPointList fConvexHullExtremePoints; ///< The points furthest from each other on hull
396  reco::ConvexHullKinkTupleList fConvexHullKinkPoints; ///< The points with large kinks along the convex hull
397 };
398 
399 /**
400  * @brief Class wrapping the above and containing volatile information to characterize the cluster
401  */
402 class ClusterParameters;
403 using ClusterParametersList = std::list<ClusterParameters>;
404 
406 {
407 public:
409  {
410  fClusterParams.clear();
411  fHitPairListPtr.clear();
412  fHit2DToHit3DListMap.clear();
413  fHit3DToEdgeMap.clear();
414  fBestHitPairListPtr.clear();
415  fBestEdgeList.clear();
416  fConvexHull.clear();
417  fFaceList.clear();
418  fVertexList.clear();
419  fHalfEdgeList.clear();
420  fClusterParameters.clear();
421  }
422 
423  ClusterParameters(reco::HitPairClusterMap::iterator& mapItr) : fHitPairListPtr(mapItr->second)
424  {
425  fClusterParams.clear();
426  fHit2DToHit3DListMap.clear();
427  fHit3DToEdgeMap.clear();
428  fBestHitPairListPtr.clear();
429  fBestEdgeList.clear();
430  fConvexHull.clear();
431  fFaceList.clear();
432  fVertexList.clear();
433  fHalfEdgeList.clear();
434  }
435 
437  {
438  fClusterParams.clear();
439  fHit2DToHit3DListMap.clear();
440  fHit3DToEdgeMap.clear();
441  fBestHitPairListPtr.clear();
442  fBestEdgeList.clear();
443  fConvexHull.clear();
444  fFaceList.clear();
445  fVertexList.clear();
446  fHalfEdgeList.clear();
447  }
448 
450 
452  {
453  fClusterParams[hit->WireID().Plane].UpdateParameters(hit);
454  }
455 
456  void addHit3D(const reco::ClusterHit3D* hit3D)
457  {
458  fHitPairListPtr.emplace_back(hit3D);
459 
460  for(const auto& hit2D : hit3D->getHits())
461  if (hit2D) fHit2DToHit3DListMap[hit2D].emplace_back(hit3D);
462  }
463 
465  {
466  for(const auto& hit3D : fHitPairListPtr)
467  {
468  for(const auto& hit2D : hit3D->getHits())
469  if (hit2D) fHit2DToHit3DListMap[hit2D].emplace_back(hit3D);
470  }
471  }
472 
485 
486  friend bool operator < (const ClusterParameters &a, const ClusterParameters& b)
487  {
488  return a.fHitPairListPtr.size() > b.fHitPairListPtr.size();
489  }
490 
491 private:
493  reco::HitPairListPtr fHitPairListPtr; // This contains the list of 3D hits in the cluster
494  reco::Hit2DToHit3DListMap fHit2DToHit3DListMap; // Provides a mapping between 2D hits and 3D hits they make
495  reco::PrincipalComponents fFullPCA; // PCA run over full set of 3D hits
496  reco::PrincipalComponents fSkeletonPCA; // PCA run over just the "skeleton" 3D hits
499  reco::EdgeList fBestEdgeList; // This has become multiuse... really need to split it up
500  reco::ConvexHull fConvexHull; // Convex hull object
501  dcel2d::FaceList fFaceList; // Keeps track of "faces" from Voronoi Diagram
502  dcel2d::VertexList fVertexList; // Keeps track of "vertices" from Voronoi Diagram
503  dcel2d::HalfEdgeList fHalfEdgeList; // Keeps track of "halfedges" from Voronoi Diagram
504  ClusterParametersList fClusterParameters; // For possible daughter clusters
505 };
506 
507 using ClusterToHitPairSetPair = std::pair<reco::ClusterParameters*,HitPairSetPtr>;
508 using ClusterToHitPairSetMap = std::unordered_map<reco::ClusterParameters*,HitPairSetPtr>;
509 using Hit2DToHit3DSetMap = std::unordered_map<const reco::ClusterHit2D*,HitPairSetPtr>;
510 using Hit2DToClusterMap = std::unordered_map<const reco::ClusterHit2D*,ClusterToHitPairSetMap>;
511 
512 }
513 
514 #endif //RECO_CLUSTER3D_H
void initialize(size_t id, unsigned int statusBits, const Eigen::Vector3f &position, float totalCharge, float avePeakTime, float deltaPeakTime, float sigmaPeakTime, float hitChiSquare, float overlapFraction, float chargeAsymmetry, float docaToAxis, float arclenToPoca, const ClusterHit2DVec &hitVec, const std::vector< float > &hitDelTSigVec, const std::vector< geo::WireID > &wireIDVec)
Definition: Cluster3D.cxx:137
Hit contains 2D hit from view 2 (w plane)
Definition: Cluster3D.h:115
reco::HitPairListPtr & getBestHitPairListPtr()
Definition: Cluster3D.h:479
void setHit(const recob::Hit *hit)
Definition: Cluster3D.h:84
std::map< int, HitPairListPtr > HitPairClusterMap
Definition: Cluster3D.h:337
PrincipalComponents m_pcaResults
Output of the prinicipal componenets analysis.
Definition: Cluster3D.h:267
reco::ConvexHullKinkTupleList fConvexHullKinkPoints
The points with large kinks along the convex hull.
Definition: Cluster3D.h:396
void setAveHitDoca(double doca) const
Definition: Cluster3D.h:251
float fDeltaPeakTime
Largest delta peak time of associated recob::Hits.
Definition: Cluster3D.h:208
void setWireID(const geo::WireID &wid) const
Definition: Cluster3D.cxx:172
reco::HitPairListPtr fBestHitPairListPtr
Definition: Cluster3D.h:498
std::list< const reco::ClusterHit2D * > Hit2DListPtr
export some data structure definitions
Definition: Cluster3D.h:333
bool m_svdOK
SVD Decomposition was successful.
Definition: Cluster3D.h:232
Cluster3D operator+(Cluster3D)
Definition: Cluster3D.cxx:268
std::list< reco::ClusterHit3D > HitPairList
Definition: Cluster3D.h:338
dcel2d::HalfEdgeList fHalfEdgeList
Definition: Cluster3D.h:503
bool getSvdOK() const
Definition: Cluster3D.h:243
unsigned m_statusBits
Default constructor.
Definition: Cluster3D.h:266
reco::PrincipalComponents fSkeletonPCA
Definition: Cluster3D.h:496
float getTotalCharge() const
Definition: Cluster3D.h:283
void clearStatusBits(unsigned bits) const
Definition: Cluster3D.h:80
float getTimeTicks() const
Definition: Cluster3D.h:75
float fTotalCharge
Sum of charges of all associated recob::Hits.
Definition: Cluster3D.h:206
void flipAxis(size_t axis)
Definition: Cluster3D.cxx:208
std::list< HalfEdge > HalfEdgeList
Definition: DCEL.h:184
double m_aveHitDoca
Average doca of hits used in PCA.
Definition: Cluster3D.h:237
enum geo::_plane_proj View_t
Enumerate the possible plane projections.
const Eigen::Vector3f getPosition() const
Definition: Cluster3D.h:157
std::list< ProjectedPoint > ProjectedPointList
Definition: Cluster3D.h:352
3D hit has 2D hits shared between clusters
Definition: Cluster3D.h:109
void setArcLenToPoca(float poca) const
Definition: Cluster3D.h:82
const bool operator<(const reco::ClusterHit3D &other) const
Definition: Cluster3D.h:187
reco::PrincipalComponents & getSkeletonPCA()
Definition: Cluster3D.h:477
Hit labelled outlier in PCA.
Definition: Cluster3D.h:112
void clearStatusBits(unsigned bits) const
Definition: Cluster3D.h:180
The data type to uniquely identify a Plane.
Definition: geo_types.h:472
void setArclenToPoca(double poca) const
Definition: Cluster3D.h:182
const geo::WireID & WireID() const
Definition: Cluster3D.h:76
void setDocaToAxis(double doca) const
Definition: Cluster3D.h:181
reco::Hit2DToHit3DListMap & getHit2DToHit3DListMap()
Definition: Cluster3D.h:474
const float * getStartPosition() const
Definition: Cluster3D.h:284
unsigned m_statusBits
Volatile status information of this 3D hit.
Definition: Cluster3D.h:42
dcel2d::VertexList fVertexList
Definition: Cluster3D.h:502
Hit is an &quot;edge&quot; hit.
Definition: Cluster3D.h:100
size_t fID
&quot;id&quot; of this hit (useful for indexing)
Definition: Cluster3D.h:203
Labelled &quot;noise&quot; by a clustering algorithm.
Definition: Cluster3D.h:107
const std::vector< float > getHitDelTSigVec() const
Definition: Cluster3D.h:171
reco::ConvexHull fConvexHull
Definition: Cluster3D.h:500
float getTotalCharge() const
Definition: Cluster3D.h:161
Hit has been rejected for any reason.
Definition: Cluster3D.h:98
reco::EdgeList & getBestEdgeList()
Definition: Cluster3D.h:480
Hit contains 2D hit from view 0 (u plane)
Definition: Cluster3D.h:113
int getNumHitsUsed() const
Definition: Cluster3D.h:244
reco::Hit3DToEdgeMap & getHit3DToEdgeMap()
Definition: Cluster3D.h:478
std::list< HitPairListPtr > HitPairListPtrList
Definition: Cluster3D.h:336
A utility class used in construction of 3D clusters.
Definition: Cluster3D.h:299
reco::HitPairListPtr & getHitPairListPtr()
Definition: Cluster3D.h:475
float getSigmaPeakTime() const
Definition: Cluster3D.h:164
process_name hit
Definition: cheaterreco.fcl:51
float fSigmaPeakTime
Quad sum of peak time sigmas.
Definition: Cluster3D.h:209
float m_docaToAxis
DOCA of hit at POCA to associated cluster axis.
Definition: Cluster3D.h:43
friend bool operator<(const ClusterHit2D &a, const ClusterHit2D &b)
Definition: Cluster3D.cxx:57
reco::PlaneToClusterParamsMap & getClusterParams()
Definition: Cluster3D.h:473
float fAvePeakTime
Average peak time of all associated recob::Hits.
Definition: Cluster3D.h:207
float getOverlapFraction() const
Definition: Cluster3D.h:166
reco::Hit3DToEdgeMap fHit3DToEdgeMap
Definition: Cluster3D.h:497
unsigned int getStatusBits() const
Definition: Cluster3D.h:156
const bool operator==(const reco::ClusterHit3D &other) const
Definition: Cluster3D.h:193
float getDocaToAxis() const
Definition: Cluster3D.h:168
friend bool operator<(const Cluster3D &a, const Cluster3D &b)
Definition: Cluster3D.cxx:342
const recob::Hit * getHit() const
Definition: Cluster3D.h:77
Hit has been used in Cluster Splitting MST.
Definition: Cluster3D.h:111
std::unordered_map< const reco::ClusterHit2D *, HitPairSetPtr > Hit2DToHit3DSetMap
Definition: Cluster3D.h:509
IDparameter< geo::PlaneID > PlaneID
Member type of validated geo::PlaneID parameter.
process_name gaushit a
float getY() const
Definition: Cluster3D.h:159
unsigned getStatusBits() const
Definition: Cluster3D.h:281
std::vector< geo::WireID > fWireIDVector
Wire ID&#39;s for the planes making up hit.
Definition: Cluster3D.h:217
ClusterParametersList & daughterList()
Definition: Cluster3D.h:449
float getAvePeakTime() const
Definition: Cluster3D.h:162
float getDocaToAxis() const
Definition: Cluster3D.h:72
const float * getEndPosition() const
Definition: Cluster3D.h:285
std::vector< float > fHitDelTSigVec
Delta t of hit to matching pair / sig.
Definition: Cluster3D.h:216
std::pair< const reco::ClusterHit3D *, reco::EdgeList > Hit3DToEdgePair
Definition: Cluster3D.h:345
std::list< EdgeTuple > EdgeList
Definition: Cluster3D.h:344
EigenValues m_eigenValues
Eigen values from SVD decomposition.
Definition: Cluster3D.h:234
Point is on primary cluster convex hull.
Definition: Cluster3D.h:103
ClusterParametersList fClusterParameters
Definition: Cluster3D.h:504
const recob::Hit * m_hit
Hit we are augmenting.
Definition: Cluster3D.h:48
const EigenValues & getEigenValues() const
Definition: Cluster3D.h:245
process_name standard_reco_uboone reco
PlaneToClusterParamsMap fClusterParams
Definition: Cluster3D.h:492
reco::PrincipalComponents & getFullPCA()
Definition: Cluster3D.h:476
reco::Hit3DToEdgeMap & getConvexHullEdgeMap()
Definition: Cluster3D.h:384
float getChargeAsymmetry() const
Definition: Cluster3D.h:167
float m_arcLenToPoca
arc length to POCA along cluster axis
Definition: Cluster3D.h:44
std::list< Face > FaceList
Definition: DCEL.h:183
float m_totalCharge
Total charge in the cluster.
Definition: Cluster3D.h:268
Eigen::Vector3f EigenValues
Definition: Cluster3D.h:225
unsigned int fStatusBits
Volatile status information of this 3D hit.
Definition: Cluster3D.h:204
reco::ProjectedPointList & getConvexHullExtremePoints()
Definition: Cluster3D.h:386
float m_timeTicks
The time (in ticks) for this hit.
Definition: Cluster3D.h:46
reco::ConvexHullKinkTupleList & getConvexHullKinkPoints()
Definition: Cluster3D.h:387
Define a container for working with the convex hull.
Definition: Cluster3D.h:359
std::unordered_map< const reco::ClusterHit2D *, reco::HitPairListPtr > Hit2DToHit3DListMap
Definition: Cluster3D.h:347
const Eigen::Vector3f & getAvePosition() const
Definition: Cluster3D.h:247
float fOverlapFraction
Hit overlap fraction start/stop of triplet.
Definition: Cluster3D.h:211
void clearStatusBits(unsigned bits) const
Definition: Cluster3D.h:289
reco::EdgeList fBestEdgeList
Definition: Cluster3D.h:499
Hit contains 2D hit from view 1 (v plane)
Definition: Cluster3D.h:114
Path checking algorithm has seen this hit.
Definition: Cluster3D.h:110
std::list< ConvexHullKinkTuple > ConvexHullKinkTupleList
Definition: Cluster3D.h:354
float m_endPosition[3]
&quot;end&quot; position for cluster
Definition: Cluster3D.h:270
size_t getID() const
Definition: Cluster3D.h:155
Eigen::Vector3f m_avePosition
Average position of hits fed to PCA.
Definition: Cluster3D.h:236
std::tuple< const reco::ClusterHit3D *, const reco::ClusterHit3D *, double > EdgeTuple
Definition: Cluster3D.h:343
std::list< const reco::ClusterHit3D * > HitPairListPtr
Definition: Cluster3D.h:334
void UpdateParameters(const reco::ClusterHit2D *hit)
Definition: Cluster3D.h:451
void setDocaToAxis(float doca) const
Definition: Cluster3D.h:81
float getXPosition() const
Definition: Cluster3D.h:74
reco::Hit2DToHit3DListMap fHit2DToHit3DListMap
Definition: Cluster3D.h:494
float fHitChiSquare
Hit ChiSquare relative to the average time.
Definition: Cluster3D.h:210
PlaneID_t Plane
Index of the plane within its TPC.
Definition: geo_types.h:493
float fDocaToAxis
DOCA to the associated cluster axis.
Definition: Cluster3D.h:213
std::pair< reco::ClusterParameters *, HitPairSetPtr > ClusterToHitPairSetPair
Definition: Cluster3D.h:507
Definition of data types for geometry description.
friend bool operator<(const ClusterParameters &a, const ClusterParameters &b)
Definition: Cluster3D.h:486
reco::EdgeList fConvexHullEdgeList
An edge list translated back to 3D hits.
Definition: Cluster3D.h:394
std::vector< const reco::ClusterHit2D * > ClusterHit2DVec
Definition: Cluster3D.h:91
std::tuple< ProjectedPoint, Eigen::Vector2f, Eigen::Vector2f > ConvexHullKinkTuple
Point plus edges that point to it.
Definition: Cluster3D.h:353
Is a convex hull extreme point.
Definition: Cluster3D.h:104
std::pair< reco::PrincipalComponents, reco::HitPairClusterMap::iterator > PCAHitPairClusterMapPair
Definition: Cluster3D.h:341
reco::ConvexHull & getConvexHull()
Definition: Cluster3D.h:481
reco::ProjectedPointList fConvexHullPointList
The points on the convex hull.
Definition: Cluster3D.h:392
reco::Hit3DToEdgeMap fConvexHullEdgeMap
Map from 3D hit to associated edge.
Definition: Cluster3D.h:393
ClusterHit2DVec & getHits()
Definition: Cluster3D.h:174
void setID(const size_t &id) const
Definition: Cluster3D.h:178
bool bitsAreSet(const unsigned int &bitsToCheck) const
Definition: Cluster3D.h:176
const PrincipalComponents & getPcaResults() const
Definition: Cluster3D.h:282
float fChargeAsymmetry
Assymetry of average of two closest to third charge.
Definition: Cluster3D.h:212
const float getAveHitDoca() const
Definition: Cluster3D.h:248
float m_startPosition[3]
&quot;start&quot; position for cluster (world coordinates)
Definition: Cluster3D.h:269
Hit has been made into Space Point.
Definition: Cluster3D.h:102
dcel2d::FaceList & getFaceList()
Definition: Cluster3D.h:482
ClusterHit2DVec fHitVector
Hits comprising this 3D hit.
Definition: Cluster3D.h:215
float fArclenToPoca
arc length along axis to DOCA point
Definition: Cluster3D.h:214
geo::WireID m_wireID
Keep track this particular hit&#39;s wireID.
Definition: Cluster3D.h:47
const std::vector< geo::WireID > & getWireIDs() const
Definition: Cluster3D.h:172
Eigen::Matrix3f EigenVectors
Definition: Cluster3D.h:226
Hit is part of Seed for track fits.
Definition: Cluster3D.h:101
float getArclenToPoca() const
Definition: Cluster3D.h:169
std::map< size_t, RecobClusterParameters > PlaneToClusterParamsMap
Definition: Cluster3D.h:342
float getDeltaPeakTime() const
Definition: Cluster3D.h:163
int m_clusterIdx
ID for this cluster.
Definition: Cluster3D.h:271
std::unordered_map< const reco::ClusterHit2D *, ClusterToHitPairSetMap > Hit2DToClusterMap
Definition: Cluster3D.h:510
reco::ProjectedPointList fProjectedPointList
The input set of points projected onto plane encompassed by the hull.
Definition: Cluster3D.h:391
Skeleton hit position averaged.
Definition: Cluster3D.h:105
ClusterParameters(reco::HitPairClusterMap::iterator &mapItr)
Definition: Cluster3D.h:423
reco::PrincipalComponents fFullPCA
Definition: Cluster3D.h:495
ClusterHit2DVec m_hitVector
Definition: Cluster3D.h:326
dcel2d::FaceList fFaceList
Definition: Cluster3D.h:501
2D representation of charge deposited in the TDC/wire plane
Definition: Hit.h:48
friend bool operator<(const PrincipalComponents &a, const PrincipalComponents &b)
Definition: Cluster3D.cxx:235
EigenVectors m_eigenVectors
The three principle axes.
Definition: Cluster3D.h:235
float getHitChiSquare() const
Definition: Cluster3D.h:165
float m_xPosition
The x coordinate for this hit.
Definition: Cluster3D.h:45
std::unordered_map< const reco::ClusterHit3D *, reco::EdgeList > Hit3DToEdgeMap
Definition: Cluster3D.h:346
&quot;visited&quot; by a clustering algorithm
Definition: Cluster3D.h:106
reco::ProjectedPointList & getConvexHullPointList()
Definition: Cluster3D.h:383
Hit is a &quot;skeleton&quot; hit.
Definition: Cluster3D.h:99
reco::ProjectedPointList & getProjectedPointList()
Definition: Cluster3D.h:382
friend std::ostream & operator<<(std::ostream &o, const ClusterHit2D &c)
Definition: Cluster3D.cxx:50
void setStatusBit(unsigned bits) const
Definition: Cluster3D.h:288
std::unordered_map< reco::ClusterParameters *, HitPairSetPtr > ClusterToHitPairSetMap
Definition: Cluster3D.h:508
reco::EdgeList & getConvexHullEdgeList()
Definition: Cluster3D.h:385
const ClusterHit2DVec & getHits() const
Definition: Cluster3D.h:170
float getZ() const
Definition: Cluster3D.h:160
std::tuple< float, float, const reco::ClusterHit3D * > ProjectedPoint
Projected coordinates and pointer to hit.
Definition: Cluster3D.h:351
ClusterParameters(reco::HitPairListPtr &hitList)
Definition: Cluster3D.h:436
void setPosition(const Eigen::Vector3f &pos) const
Definition: Cluster3D.h:185
dcel2d::HalfEdgeList & getHalfEdgeList()
Definition: Cluster3D.h:484
std::list< Vertex > VertexList
Definition: DCEL.h:182
reco::HitPairListPtr fHitPairListPtr
Definition: Cluster3D.h:493
int getClusterIdx() const
Definition: Cluster3D.h:286
std::list< ClusterParameters > ClusterParametersList
Definition: Cluster3D.h:403
void UpdateParameters(const reco::ClusterHit2D *hit)
Definition: Cluster3D.cxx:361
const EigenVectors & getEigenVectors() const
Definition: Cluster3D.h:246
attached to a cluster
Definition: Cluster3D.h:108
friend std::ostream & operator<<(std::ostream &o, const Cluster3D &c)
Definition: Cluster3D.cxx:322
void addHit3D(const reco::ClusterHit3D *hit3D)
Definition: Cluster3D.h:456
dcel2d::VertexList & getVertexList()
Definition: Cluster3D.h:483
std::set< const reco::ClusterHit3D * > HitPairSetPtr
Definition: Cluster3D.h:335
Eigen::Vector3f fPosition
position of this hit combination in world coordinates
Definition: Cluster3D.h:205
friend std::ostream & operator<<(std::ostream &o, const ClusterHit3D &c)
Definition: Cluster3D.cxx:177
float getArcLenToPoca() const
Definition: Cluster3D.h:73
float getX() const
Definition: Cluster3D.h:158
void fillHit2DToHit3DListMap()
Definition: Cluster3D.h:464
reco::ProjectedPointList fConvexHullExtremePoints
The points furthest from each other on hull.
Definition: Cluster3D.h:395
unsigned getStatusBits() const
Definition: Cluster3D.h:71
void setStatusBit(unsigned bits) const
Definition: Cluster3D.h:179
void setStatusBit(unsigned bits) const
Definition: Cluster3D.h:79
int m_numHitsUsed
Number of hits in the decomposition.
Definition: Cluster3D.h:233
friend std::ostream & operator<<(std::ostream &o, const PrincipalComponents &a)
Definition: Cluster3D.cxx:215