LArSoft  v06_85_00
Liquid Argon Software toolkit - http://larsoft.org/
Cluster3D.h
Go to the documentation of this file.
1 
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 #include <memory>
24 
28 
29 // Eigen
30 #include <Eigen/Dense>
31 
32 namespace reco {
33 
34 // First define a container object to augment the sparse 2D hit information
36 {
37 public:
38 
39  ClusterHit2D(); // Default constructor
40 
41 private:
42 
43  mutable unsigned m_statusBits;
44  mutable float m_docaToAxis;
45  mutable float m_arcLenToPoca;
46  float m_xPosition;
47  float m_timeTicks;
48  const recob::Hit* m_hit;
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 recob::Hit& recobHit);
67 
68  unsigned getStatusBits() const {return m_statusBits;}
69  float getDocaToAxis() const {return m_docaToAxis;}
70  float getArcLenToPoca() const {return m_arcLenToPoca;}
71  float getXPosition() const {return m_xPosition;}
72  float getTimeTicks() const {return m_timeTicks;}
73  const recob::Hit& getHit() const {return *m_hit;}
74 
75  void setStatusBit(unsigned bits) const {m_statusBits |= bits;}
76  void clearStatusBits(unsigned bits) const {m_statusBits &= ~bits;}
77  void setDocaToAxis(float doca) const {m_docaToAxis = doca;}
78  void setArcLenToPoca(float poca) const {m_arcLenToPoca = poca;}
79 
80  friend std::ostream& operator << (std::ostream& o, const ClusterHit2D& c);
81  friend bool operator < (const ClusterHit2D & a, const ClusterHit2D & b);
82 
83 };
84 
85 using ClusterHit2DVec = std::vector<const reco::ClusterHit2D*>;
86 
87 // Now define an object with the recob::Hit information that will comprise the 3D cluster
89 {
90 public:
91 
92  enum StatusBits { REJECTEDHIT = 0x80000000,
93  SKELETONHIT = 0x10000000,
94  EDGEHIT = 0x20000000,
95  SEEDHIT = 0x40000000,
96  MADESPACEPOINT = 0x08000000,
97  SKELETONPOSAVE = 0x00100000,
98  CLUSTERVISITED = 0x00008000,
99  CLUSTERNOISE = 0x00004000,
100  CLUSTERATTACHED = 0x00002000,
101  CLUSTERSHARED = 0x00001000,
102  PATHCHECKED = 0x00000800,
103  SELECTEDBYMST = 0x00000100,
104  PCAOUTLIER = 0x00000080,
105  HITINVIEW0 = 0x00000001,
106  HITINVIEW1 = 0x00000002,
107  HITINVIEW2 = 0x00000004
108  };
109 
110 
111  ClusterHit3D(); // Default constructor
112 
113  ClusterHit3D(size_t id,
114  unsigned int statusBits,
115  const float* position,
116  float totalCharge,
117  float avePeakTime,
118  float deltaPeakTime,
119  float sigmaPeakTime,
120  float hitChiSquare,
121  float docaToAxis,
122  float arclenToPoca,
123  const ClusterHit2DVec& hitVec,
124  const std::vector<float>& hitDelTSigVec,
125  const std::vector<geo::WireID>& wireIDVec);
126 
127  ClusterHit3D(const ClusterHit3D&);
128 
129  void initialize(size_t id,
130  unsigned int statusBits,
131  const float* position,
132  float totalCharge,
133  float avePeakTime,
134  float deltaPeakTime,
135  float sigmaPeakTime,
136  float hitChiSquare,
137  float docaToAxis,
138  float arclenToPoca,
139  const ClusterHit2DVec& hitVec,
140  const std::vector<float>& hitDelTSigVec,
141  const std::vector<geo::WireID>& wireIDVec);
142 
143  size_t getID() const {return m_id;}
144  unsigned int getStatusBits() const {return m_statusBits;}
145  const float* getPosition() const {return m_position;}
146  float getX() const {return m_position[0];}
147  float getY() const {return m_position[1];}
148  float getZ() const {return m_position[2];}
149  float getTotalCharge() const {return m_totalCharge;}
150  float getAvePeakTime() const {return m_avePeakTime;}
151  float getDeltaPeakTime() const {return m_deltaPeakTime;}
152  float getSigmaPeakTime() const {return m_sigmaPeakTime;}
153  float getHitChiSquare() const {return m_hitChiSquare;}
154  float getDocaToAxis() const {return m_docaToAxis;}
155  float getArclenToPoca() const {return m_arclenToPoca;}
156  const ClusterHit2DVec& getHits() const {return m_hitVector;}
157  const std::vector<float> getHitDelTSigVec() const {return m_hitDelTSigVec;}
158  const std::vector<geo::WireID>& getWireIDs() const {return m_wireIDVector;}
159 
160  bool bitsAreSet(const unsigned int& bitsToCheck) const {return m_statusBits & bitsToCheck;}
161 
162  void setID(const size_t& id) const {m_id = id;}
163  void setStatusBit(unsigned bits) const {m_statusBits |= bits;}
164  void clearStatusBits(unsigned bits) const {m_statusBits &= ~bits;}
165  void setDocaToAxis(double doca) const {m_docaToAxis = doca;}
166  void setArclenToPoca(double poca) const {m_arclenToPoca = poca;}
167  void setWireID(const geo::WireID& wid) const;
168 
169  void setPosition(const float* pos) const {m_position[0] = pos[0]; m_position[1] = pos[1]; m_position[2] = pos[2];}
170 
171  const bool operator<(const reco::ClusterHit3D& other) const
172  {
173  if (m_position[2] != other.m_position[2]) return m_position[2] < other.m_position[2];
174  else return m_position[0] < other.m_position[0];
175  }
176 
177  const bool operator==(const reco::ClusterHit3D& other) const
178  {
179  return m_id == other.m_id;
180  }
181 
182  friend std::ostream& operator << (std::ostream& o, const ClusterHit3D& c);
183  //friend bool operator < (const ClusterHit3D & a, const ClusterHit3D & b);
184 
185 private:
186 
187  mutable size_t m_id;
188  mutable unsigned int m_statusBits;
189  mutable float m_position[3];
195  mutable float m_docaToAxis;
196  mutable float m_arclenToPoca;
198  mutable std::vector<float> m_hitDelTSigVec;
199  mutable std::vector<geo::WireID> m_wireIDVector;
200 };
201 
202 // We also need to define a container for the output of the PCA Analysis
204 {
205 public:
206 
207  typedef std::vector<std::vector<float>> EigenVectors;
208 
210 
211 private:
212 
213  bool m_svdOK;
215  float m_eigenValues[3];
216  EigenVectors m_eigenVectors;
217  float m_avePosition[3];
218  mutable double m_aveHitDoca;
219 
220 public:
221 
222  PrincipalComponents(bool ok, int nHits, const float* eigenValues, const EigenVectors& eigenVecs, const float* avePos, const float aveHitDoca = 9999.);
223 
224  bool getSvdOK() const {return m_svdOK;}
225  int getNumHitsUsed() const {return m_numHitsUsed;}
226  const float* getEigenValues() const {return m_eigenValues;}
227  const EigenVectors& getEigenVectors() const {return m_eigenVectors;}
228  const float* getAvePosition() const {return m_avePosition;}
229  const float getAveHitDoca() const {return m_aveHitDoca;}
230 
231  void flipAxis(size_t axis);
232  void setAveHitDoca(double doca) const {m_aveHitDoca = doca;}
233 
234  friend std::ostream& operator << (std::ostream & o, const PrincipalComponents& a);
235  friend bool operator < (const PrincipalComponents& a, const PrincipalComponents& b);
236 
237 };
238 
240 {
241 public:
242 
243  Cluster3D();
244 
245 private:
246 
247  mutable unsigned m_statusBits;
250  float m_startPosition[3];
251  float m_endPosition[3];
253 
254 public:
255  Cluster3D(unsigned statusBits,
256  const PrincipalComponents& pcaResults,
257  float totalCharge,
258  const float* startPosition,
259  const float* endPosition,
260  int idx);
261 
262  unsigned getStatusBits() const {return m_statusBits;}
263  const PrincipalComponents& getPcaResults() const {return m_pcaResults;}
264  float getTotalCharge() const {return m_totalCharge;}
265  const float* getStartPosition() const {return m_startPosition;}
266  const float* getEndPosition() const {return m_endPosition;}
267  int getClusterIdx() const {return m_clusterIdx;}
268 
269  void setStatusBit(unsigned bits) const {m_statusBits |= bits;}
270  void clearStatusBits(unsigned bits) const {m_statusBits &= ~bits;}
271 
273  friend std::ostream& operator << (std::ostream& o, const Cluster3D& c);
274  friend bool operator < (const Cluster3D & a, const Cluster3D & b);
275 };
276 
281 {
282 public:
283  RecobClusterParameters() : m_startTime(999999.),
284  m_sigmaStartTime(1.),
285  m_endTime(0.),
286  m_sigmaEndTime(1.),
287  m_totalCharge(0.),
288  m_startWire(9999999),
289  m_endWire(0),
290  m_plane(100),
291  m_view(geo::kUnknown)
292  {
293  m_hitVector.clear();
294  }
295 
296  void UpdateParameters(const reco::ClusterHit2D* hit);
297 
298  float m_startTime;
300  float m_endTime;
303  unsigned int m_startWire;
304  unsigned int m_endWire;
305  unsigned int m_plane;
308 };
309 
310 
314 using Hit2DListPtr = std::list<const reco::ClusterHit2D*>;
315 using HitPairListPtr = std::list<const reco::ClusterHit3D*>;
316 using HitPairSetPtr = std::set<const reco::ClusterHit3D*>;
317 using HitPairListPtrList = std::list<HitPairListPtr>;
318 using HitPairClusterMap = std::map<int, HitPairListPtr>;
319 using HitPairList = std::list<std::unique_ptr<reco::ClusterHit3D>>;
320 
321 using PCAHitPairClusterMapPair = std::pair<reco::PrincipalComponents, reco::HitPairClusterMap::iterator>;
322 using PlaneToClusterParamsMap = std::map<size_t, RecobClusterParameters>;
323 using EdgeTuple = std::tuple<const reco::ClusterHit3D*,const reco::ClusterHit3D*,double>;
324 using EdgeList = std::list<EdgeTuple>;
325 using Hit3DToEdgePair = std::pair<const reco::ClusterHit3D*, reco::EdgeList>;
326 using Hit3DToEdgeMap = std::unordered_map<const reco::ClusterHit3D*, reco::EdgeList>;
327 using Hit2DToHit3DListMap = std::unordered_map<const reco::ClusterHit2D*, reco::HitPairListPtr>;
328 //using VertexPoint = Eigen::Vector3f;
329 //using VertexPointList = std::list<Eigen::Vector3f>;
330 
334 class ClusterParameters;
335 using ClusterParametersList = std::list<ClusterParameters>;
336 
338 {
339 public:
341  {
342  m_clusterParams.clear();
343  m_hitPairListPtr.clear();
344  m_hit2DToHit3DListMap.clear();
345  m_hit3DToEdgeMap.clear();
346  m_bestHitPairListPtr.clear();
347  m_bestEdgeList.clear();
348  m_faceList.clear();
349  m_vertexList.clear();
350  m_halfEdgeList.clear();
351  m_clusterParameters.clear();
352  }
353 
354  ClusterParameters(reco::HitPairClusterMap::iterator& mapItr) : m_hitPairListPtr(mapItr->second)
355  {
356  m_clusterParams.clear();
357  m_hit2DToHit3DListMap.clear();
358  m_hit3DToEdgeMap.clear();
359  m_bestHitPairListPtr.clear();
360  m_bestEdgeList.clear();
361  m_faceList.clear();
362  m_vertexList.clear();
363  m_halfEdgeList.clear();
364  }
365 
366  ClusterParameters(reco::HitPairListPtr& hitList) : m_hitPairListPtr(hitList)
367  {
368  m_clusterParams.clear();
369  m_hit2DToHit3DListMap.clear();
370  m_hit3DToEdgeMap.clear();
371  m_bestHitPairListPtr.clear();
372  m_bestEdgeList.clear();
373  m_faceList.clear();
374  m_vertexList.clear();
375  m_halfEdgeList.clear();
376  }
377 
378  ClusterParametersList& daughterList() {return m_clusterParameters;}
379 
381  {
382  m_clusterParams[hit->getHit().WireID().Plane].UpdateParameters(hit);
383  }
384 
385  void addHit3D(const reco::ClusterHit3D* hit3D)
386  {
387  m_hitPairListPtr.emplace_back(hit3D);
388 
389  for(const auto& hit2D : hit3D->getHits())
390  if (hit2D) m_hit2DToHit3DListMap[hit2D].emplace_back(hit3D);
391  }
392 
394  {
395  for(const auto& hit3D : m_hitPairListPtr)
396  {
397  for(const auto& hit2D : hit3D->getHits())
398  if (hit2D) m_hit2DToHit3DListMap[hit2D].emplace_back(hit3D);
399  }
400  }
401 
402  reco::PlaneToClusterParamsMap& getClusterParams() {return m_clusterParams;}
403  reco::Hit2DToHit3DListMap& getHit2DToHit3DListMap() {return m_hit2DToHit3DListMap;}
404  reco::HitPairListPtr& getHitPairListPtr() {return m_hitPairListPtr;}
405  reco::PrincipalComponents& getFullPCA() {return m_fullPCA;}
406  reco::PrincipalComponents& getSkeletonPCA() {return m_skeletonPCA;}
407  reco::Hit3DToEdgeMap& getHit3DToEdgeMap() {return m_hit3DToEdgeMap;}
408  reco::HitPairListPtr& getBestHitPairListPtr() {return m_bestHitPairListPtr;}
409  reco::EdgeList& getBestEdgeList() {return m_bestEdgeList;}
410  dcel2d::FaceList& getFaceList() {return m_faceList;}
411  dcel2d::VertexList& getVertexList() {return m_vertexList;}
412  dcel2d::HalfEdgeList& getHalfEdgeList() {return m_halfEdgeList;}
413 
414  friend bool operator < (const ClusterParameters &a, const ClusterParameters& b)
415  {
416  return a.m_hitPairListPtr.size() > b.m_hitPairListPtr.size();
417  }
418 
419 private:
421  reco::HitPairListPtr m_hitPairListPtr; // This contains the list of 3D hits in the cluster
422  reco::Hit2DToHit3DListMap m_hit2DToHit3DListMap; // Provides a mapping between 2D hits and 3D hits they make
423  reco::PrincipalComponents m_fullPCA; // PCA run over full set of 3D hits
424  reco::PrincipalComponents m_skeletonPCA; // PCA run over just the "skeleton" 3D hits
428  dcel2d::FaceList m_faceList; // Keeps track of "faces" from Voronoi Diagram
429  dcel2d::VertexList m_vertexList; // Keeps track of "vertices" from Voronoi Diagram
430  dcel2d::HalfEdgeList m_halfEdgeList; // Keeps track of "halfedges" from Voronoi Diagram
431  ClusterParametersList m_clusterParameters; // For possible daughter clusters
432 };
433 
434 using ClusterToHitPairSetPair = std::pair<reco::ClusterParameters*,HitPairSetPtr>;
435 using ClusterToHitPairSetMap = std::unordered_map<reco::ClusterParameters*,HitPairSetPtr>;
436 using Hit2DToHit3DSetMap = std::unordered_map<const reco::ClusterHit2D*,HitPairSetPtr>;
437 using Hit2DToClusterMap = std::unordered_map<const reco::ClusterHit2D*,ClusterToHitPairSetMap>;
438 
439 }
440 
441 #endif //RECO_CLUSTER3D_H
reco::HitPairListPtr & getBestHitPairListPtr()
Definition: Cluster3D.h:408
dcel2d::VertexList m_vertexList
Definition: Cluster3D.h:429
std::map< int, HitPairListPtr > HitPairClusterMap
Definition: Cluster3D.h:318
PrincipalComponents m_pcaResults
Output of the prinicipal componenets analysis.
Definition: Cluster3D.h:248
void setAveHitDoca(double doca) const
Definition: Cluster3D.h:232
std::list< const reco::ClusterHit2D * > Hit2DListPtr
export some data structure definitions
Definition: Cluster3D.h:314
bool m_svdOK
SVD Decomposition was successful.
Definition: Cluster3D.h:213
size_t m_id
"id" of this hit (useful for indexing)
Definition: Cluster3D.h:187
bool getSvdOK() const
Definition: Cluster3D.h:224
unsigned m_statusBits
Default constructor.
Definition: Cluster3D.h:247
float getTotalCharge() const
Definition: Cluster3D.h:264
void clearStatusBits(unsigned bits) const
Definition: Cluster3D.h:76
float getTimeTicks() const
Definition: Cluster3D.h:72
float m_hitChiSquare
Hit ChiSquare relative to the average time.
Definition: Cluster3D.h:194
float m_arclenToPoca
arc length along axis to DOCA point
Definition: Cluster3D.h:196
intermediate_table::iterator iterator
std::list< HalfEdge > HalfEdgeList
Definition: DCEL.h:180
double m_aveHitDoca
Average doca of hits used in PCA.
Definition: Cluster3D.h:218
enum geo::_plane_proj View_t
Enumerate the possible plane projections.
Unknown view.
Definition: geo_types.h:83
reco::Hit2DToHit3DListMap m_hit2DToHit3DListMap
Definition: Cluster3D.h:422
geo::WireID WireID() const
Initial tdc tick for hit.
Definition: Hit.h:234
void setArcLenToPoca(float poca) const
Definition: Cluster3D.h:78
const bool operator<(const reco::ClusterHit3D &other) const
Definition: Cluster3D.h:171
reco::EdgeList m_bestEdgeList
Definition: Cluster3D.h:427
Declaration of signal hit object.
float m_deltaPeakTime
Largest delta peak time of associated recob::Hits.
Definition: Cluster3D.h:192
reco::PrincipalComponents & getSkeletonPCA()
Definition: Cluster3D.h:406
void clearStatusBits(unsigned bits) const
Definition: Cluster3D.h:164
void setArclenToPoca(double poca) const
Definition: Cluster3D.h:166
void setDocaToAxis(double doca) const
Definition: Cluster3D.h:165
dcel2d::FaceList m_faceList
Definition: Cluster3D.h:428
reco::Hit2DToHit3DListMap & getHit2DToHit3DListMap()
Definition: Cluster3D.h:403
const float * getStartPosition() const
Definition: Cluster3D.h:265
unsigned m_statusBits
Volatile status information of this 3D hit.
Definition: Cluster3D.h:43
const std::vector< float > getHitDelTSigVec() const
Definition: Cluster3D.h:157
float getTotalCharge() const
Definition: Cluster3D.h:149
reco::EdgeList & getBestEdgeList()
Definition: Cluster3D.h:409
const float * getAvePosition() const
Definition: Cluster3D.h:228
int getNumHitsUsed() const
Definition: Cluster3D.h:225
reco::Hit3DToEdgeMap & getHit3DToEdgeMap()
Definition: Cluster3D.h:407
std::list< HitPairListPtr > HitPairListPtrList
Definition: Cluster3D.h:317
A utility class used in construction of 3D clusters.
Definition: Cluster3D.h:280
reco::HitPairListPtr & getHitPairListPtr()
Definition: Cluster3D.h:404
float getSigmaPeakTime() const
Definition: Cluster3D.h:152
reco::PrincipalComponents m_skeletonPCA
Definition: Cluster3D.h:424
float m_docaToAxis
DOCA of hit at POCA to associated cluster axis.
Definition: Cluster3D.h:44
friend bool operator<(const ClusterHit2D &a, const ClusterHit2D &b)
Definition: Cluster3D.cxx:47
reco::PlaneToClusterParamsMap & getClusterParams()
Definition: Cluster3D.h:402
dcel2d::HalfEdgeList m_halfEdgeList
Definition: Cluster3D.h:430
unsigned int getStatusBits() const
Definition: Cluster3D.h:144
const bool operator==(const reco::ClusterHit3D &other) const
Definition: Cluster3D.h:177
float getDocaToAxis() const
Definition: Cluster3D.h:154
std::unordered_map< const reco::ClusterHit2D *, HitPairSetPtr > Hit2DToHit3DSetMap
Definition: Cluster3D.h:436
std::list< std::unique_ptr< reco::ClusterHit3D >> HitPairList
Definition: Cluster3D.h:319
float getY() const
Definition: Cluster3D.h:147
unsigned getStatusBits() const
Definition: Cluster3D.h:262
ClusterParametersList & daughterList()
Definition: Cluster3D.h:378
float getAvePeakTime() const
Definition: Cluster3D.h:150
float m_sigmaPeakTime
Quad sum of peak time sigmas.
Definition: Cluster3D.h:193
float getDocaToAxis() const
Definition: Cluster3D.h:69
const float * getEndPosition() const
Definition: Cluster3D.h:266
reco::HitPairListPtr m_bestHitPairListPtr
Definition: Cluster3D.h:426
std::pair< const reco::ClusterHit3D *, reco::EdgeList > Hit3DToEdgePair
Definition: Cluster3D.h:325
std::list< EdgeTuple > EdgeList
Definition: Cluster3D.h:324
const recob::Hit * m_hit
Hit we are augmenting.
Definition: Cluster3D.h:48
float m_avePeakTime
Average peak time of all associated recob::Hits.
Definition: Cluster3D.h:191
reco::PrincipalComponents & getFullPCA()
Definition: Cluster3D.h:405
float m_arcLenToPoca
arc length to POCA along cluster axis
Definition: Cluster3D.h:45
std::list< Face > FaceList
Definition: DCEL.h:179
float m_totalCharge
Total charge in the cluster.
Definition: Cluster3D.h:249
ClusterParametersList m_clusterParameters
Definition: Cluster3D.h:431
float m_docaToAxis
DOCA to the associated cluster axis.
Definition: Cluster3D.h:195
float m_timeTicks
The time (in ticks) for this hit.
Definition: Cluster3D.h:47
std::unordered_map< const reco::ClusterHit2D *, reco::HitPairListPtr > Hit2DToHit3DListMap
Definition: Cluster3D.h:327
void clearStatusBits(unsigned bits) const
Definition: Cluster3D.h:270
std::vector< geo::WireID > m_wireIDVector
Wire ID&#39;s for the planes making up hit.
Definition: Cluster3D.h:199
const recob::Hit & getHit() const
Definition: Cluster3D.h:73
size_t getID() const
Definition: Cluster3D.h:143
std::tuple< const reco::ClusterHit3D *, const reco::ClusterHit3D *, double > EdgeTuple
Definition: Cluster3D.h:323
std::list< const reco::ClusterHit3D * > HitPairListPtr
Definition: Cluster3D.h:315
void UpdateParameters(const reco::ClusterHit2D *hit)
Definition: Cluster3D.h:380
float m_totalCharge
Sum of charges of all associated recob::Hits.
Definition: Cluster3D.h:190
void setDocaToAxis(float doca) const
Definition: Cluster3D.h:77
float getXPosition() const
Definition: Cluster3D.h:71
std::vector< std::vector< float > > EigenVectors
Definition: Cluster3D.h:207
PlaneToClusterParamsMap m_clusterParams
Definition: Cluster3D.h:420
PlaneID_t Plane
Index of the plane within its TPC.
Definition: geo_types.h:258
std::pair< reco::ClusterParameters *, HitPairSetPtr > ClusterToHitPairSetPair
Definition: Cluster3D.h:434
Definition of data types for geometry description.
const float * getPosition() const
Definition: Cluster3D.h:145
Detector simulation of raw signals on wires.
std::vector< const reco::ClusterHit2D * > ClusterHit2DVec
Definition: Cluster3D.h:85
unsigned int m_statusBits
Volatile status information of this 3D hit.
Definition: Cluster3D.h:188
std::pair< reco::PrincipalComponents, reco::HitPairClusterMap::iterator > PCAHitPairClusterMapPair
Definition: Cluster3D.h:321
void setID(const size_t &id) const
Definition: Cluster3D.h:162
bool bitsAreSet(const unsigned int &bitsToCheck) const
Definition: Cluster3D.h:160
const float * getEigenValues() const
Definition: Cluster3D.h:226
const PrincipalComponents & getPcaResults() const
Definition: Cluster3D.h:263
void setPosition(const float *pos) const
Definition: Cluster3D.h:169
const float getAveHitDoca() const
Definition: Cluster3D.h:229
QuadExpr operator+(double v, const QuadExpr &e)
Definition: QuadExpr.h:37
dcel2d::FaceList & getFaceList()
Definition: Cluster3D.h:410
reco::Hit3DToEdgeMap m_hit3DToEdgeMap
Definition: Cluster3D.h:425
const std::vector< geo::WireID > & getWireIDs() const
Definition: Cluster3D.h:158
Definition of utility objects for use in the 3D clustering for LArSoft.
Definition: Cluster3D.cxx:17
float getArclenToPoca() const
Definition: Cluster3D.h:155
std::map< size_t, RecobClusterParameters > PlaneToClusterParamsMap
Definition: Cluster3D.h:322
float getDeltaPeakTime() const
Definition: Cluster3D.h:151
reco::PrincipalComponents m_fullPCA
Definition: Cluster3D.h:423
int m_clusterIdx
ID for this cluster.
Definition: Cluster3D.h:252
std::unordered_map< const reco::ClusterHit2D *, ClusterToHitPairSetMap > Hit2DToClusterMap
Definition: Cluster3D.h:437
std::vector< float > m_hitDelTSigVec
Delta t of hit to matching pair / sig.
Definition: Cluster3D.h:198
ClusterParameters(reco::HitPairClusterMap::iterator &mapItr)
Definition: Cluster3D.h:354
ClusterHit2DVec m_hitVector
Definition: Cluster3D.h:307
2D representation of charge deposited in the TDC/wire plane
Definition: Hit.h:49
EigenVectors m_eigenVectors
The three principle axes.
Definition: Cluster3D.h:216
float getHitChiSquare() const
Definition: Cluster3D.h:153
float m_xPosition
The x coordinate for this hit.
Definition: Cluster3D.h:46
std::unordered_map< const reco::ClusterHit3D *, reco::EdgeList > Hit3DToEdgeMap
Definition: Cluster3D.h:326
reco::HitPairListPtr m_hitPairListPtr
Definition: Cluster3D.h:421
friend std::ostream & operator<<(std::ostream &o, const ClusterHit2D &c)
Definition: Cluster3D.cxx:40
void setStatusBit(unsigned bits) const
Definition: Cluster3D.h:269
std::unordered_map< reco::ClusterParameters *, HitPairSetPtr > ClusterToHitPairSetMap
Definition: Cluster3D.h:435
Namespace collecting geometry-related classes utilities.
float m_position[3]
position of this hit combination in world coordinates
Definition: Cluster3D.h:189
const ClusterHit2DVec & getHits() const
Definition: Cluster3D.h:156
float getZ() const
Definition: Cluster3D.h:148
ClusterHit2DVec m_hitVector
Hits comprising this 3D hit.
Definition: Cluster3D.h:197
ClusterParameters(reco::HitPairListPtr &hitList)
Definition: Cluster3D.h:366
dcel2d::HalfEdgeList & getHalfEdgeList()
Definition: Cluster3D.h:412
std::list< Vertex > VertexList
Definition: DCEL.h:178
int getClusterIdx() const
Definition: Cluster3D.h:267
std::list< ClusterParameters > ClusterParametersList
Definition: Cluster3D.h:335
const EigenVectors & getEigenVectors() const
Definition: Cluster3D.h:227
void addHit3D(const reco::ClusterHit3D *hit3D)
Definition: Cluster3D.h:385
dcel2d::VertexList & getVertexList()
Definition: Cluster3D.h:411
std::set< const reco::ClusterHit3D * > HitPairSetPtr
Definition: Cluster3D.h:316
float getArcLenToPoca() const
Definition: Cluster3D.h:70
float getX() const
Definition: Cluster3D.h:146
void fillHit2DToHit3DListMap()
Definition: Cluster3D.h:393
unsigned getStatusBits() const
Definition: Cluster3D.h:68
void setStatusBit(unsigned bits) const
Definition: Cluster3D.h:163
void setStatusBit(unsigned bits) const
Definition: Cluster3D.h:75
int m_numHitsUsed
Number of hits in the decomposition.
Definition: Cluster3D.h:214