CMS 3D CMS Logo

PFRecHit.h
Go to the documentation of this file.
1 #ifndef ParticleFlowReco_PFRecHit_h
2 #define ParticleFlowReco_PFRecHit_h
3 
5 #include <vector>
6 #include <map>
7 #include <memory>
8 #include <iostream>
9 
12 #include "Math/GenVector/PositionVector3D.h"
15 
18 
20 
21 
22 namespace reco {
23 
32  class PFRecHit {
33 
34  public:
40 
41  struct Neighbours {
42  using Pointer = unsigned int const *;
44  Neighbours(Pointer ib, unsigned int n) : b(ib), e(ib+n){}
46  Pointer begin() const {return b;}
47  Pointer end() const {return e;}
48  unsigned int size() const { return e-b;}
49  };
50 
51  enum {
52  NONE=0
53  };
56 
57  PFRecHit(std::shared_ptr<const CaloCellGeometry> caloCell,
58  unsigned int detId, PFLayer::Layer layer,
59  float energy) : caloCell_(std::move(caloCell)), detId_(detId),
60  layer_(layer), energy_(energy){}
61 
62 
63 
65  PFRecHit(const PFRecHit& other) = default;
66  PFRecHit(PFRecHit&& other) = default;
67  PFRecHit & operator=(const PFRecHit& other) = default;
68  PFRecHit & operator=(PFRecHit&& other) = default;
69 
70 
72  ~PFRecHit()=default;
73 
74  void setEnergy( float energy) { energy_ = energy; }
75 
76 
77  void addNeighbour(short x,short y, short z, unsigned int);
78  unsigned int getNeighbour(short x,short y, short z) const;
79  void setTime( double time) { time_ = time; }
80  void setDepth( int depth) { depth_ = depth; }
81  void clearNeighbours() {
82  neighbours_.clear();
83  neighbourInfos_.clear();
85  }
86 
89  }
92  }
93 
95  return buildNeighbours(neighbours_.size());
96  }
97 
98  const std::vector<unsigned short>& neighbourInfos() {
99  return neighbourInfos_;
100  }
101 
102 
104  CaloCellGeometry const & caloCell() const { return *(caloCell_.get()); }
105  bool hasCaloCell() const { return (caloCell_ != nullptr); }
106 
108  unsigned detId() const {return detId_;}
109 
111  PFLayer::Layer layer() const { return layer_; }
112 
114  float energy() const { return energy_; }
115 
116 
118  float time() const { return time_; }
119 
121  int depth() const { return depth_; }
122 
124  double pt2() const { return energy_ * energy_ *
125  ( position().perp2()/ position().mag2());}
126 
127 
129  PositionType const& position() const { return caloCell().getPosition().basicVector(); }
130 
131  RhoEtaPhi const& positionREP() const { return caloCell().repPos(); }
132 
134  CornersVec const& getCornersXYZ() const { return caloCell().getCorners(); }
135 
136  RepCorners const& getCornersREP() const { return caloCell().getCornersREP();}
137 
138 
140  bool operator>=(const PFRecHit& rhs) const { return (energy_>=rhs.energy_); }
141 
143  bool operator> (const PFRecHit& rhs) const { return (energy_> rhs.energy_); }
144 
146  bool operator<=(const PFRecHit& rhs) const { return (energy_<=rhs.energy_); }
147 
149  bool operator< (const PFRecHit& rhs) const { return (energy_< rhs.energy_); }
150 
151 
152  private:
153 
154  Neighbours buildNeighbours(unsigned int n) const { return Neighbours(&neighbours_.front(),n);}
155 
157  std::shared_ptr<const CaloCellGeometry> caloCell_=nullptr;
158 
160  unsigned int detId_=0;
161 
164 
166  float energy_=0;
167 
169  float time_=-1;
170 
172  int depth_=0;
173 
174 
176  std::vector< unsigned int > neighbours_;
177  std::vector< unsigned short > neighbourInfos_;
178 
179  //Caching the neighbours4/8 per request of Lindsey
180  unsigned int neighbours4_ = 0;
181  unsigned int neighbours8_ = 0;
182  };
183 
184 }
185 std::ostream& operator<<(std::ostream& out, const reco::PFRecHit& hit);
186 
187 #endif
float energy_
rechit energy
Definition: PFRecHit.h:166
std::vector< unsigned int > neighbours_
indices to existing neighbours (1 common side)
Definition: PFRecHit.h:176
EZArrayFL< GlobalPoint > CornersVec
void setDepth(int depth)
Definition: PFRecHit.h:80
BaseClass::BasicVectorType BasicVectorType
Definition: Point3DBase.h:19
Pointer begin() const
Definition: PFRecHit.h:46
double pt2() const
rechit momentum transverse to the beam, squared.
Definition: PFRecHit.h:124
unsigned int neighbours4_
Definition: PFRecHit.h:180
float time() const
timing for cleaned hits
Definition: PFRecHit.h:118
unsigned detId() const
rechit detId
Definition: PFRecHit.h:108
RepCorners const & getCornersREP() const
Definition: PFRecHit.h:136
PFRecHit & operator=(const PFRecHit &other)=default
T perp2() const
Squared magnitude of transverse component.
virtual const GlobalPoint & getPosition() const
Returns the position of reference for this cell.
int depth_
depth
Definition: PFRecHit.h:172
void addNeighbour(short x, short y, short z, unsigned int)
Definition: PFRecHit.cc:5
unsigned int detId_
cell detid
Definition: PFRecHit.h:160
PFLayer::Layer layer_
rechit layer
Definition: PFRecHit.h:163
std::vector< unsigned short > neighbourInfos_
Definition: PFRecHit.h:177
std::shared_ptr< const CaloCellGeometry > caloCell_
cell geometry
Definition: PFRecHit.h:157
Pointer end() const
Definition: PFRecHit.h:47
bool operator<=(const PFRecHit &rhs) const
comparison <= operator
Definition: PFRecHit.h:146
unsigned int getNeighbour(short x, short y, short z) const
Definition: PFRecHit.cc:49
bool operator<(const PFRecHit &rhs) const
comparison < operator
Definition: PFRecHit.h:149
PositionType const & position() const
rechit cell centre x, y, z
Definition: PFRecHit.h:129
RhoEtaPhi const & repPos() const
unsigned int neighbours8_
Definition: PFRecHit.h:181
int depth() const
depth for segemntation
Definition: PFRecHit.h:121
PFLayer::Layer layer() const
rechit layer
Definition: PFRecHit.h:111
RepCorners const & getCornersREP() const
void clearNeighbours()
Definition: PFRecHit.h:81
std::ostream & operator<<(std::ostream &, BeamSpot beam)
Definition: BeamSpot.cc:71
Particle flow rechit (rechit + geometry and topology information). See clustering algorithm in PFClus...
Definition: PFRecHit.h:32
Neighbours buildNeighbours(unsigned int n) const
Definition: PFRecHit.h:154
float time_
time
Definition: PFRecHit.h:169
bool hasCaloCell() const
Definition: PFRecHit.h:105
Neighbours neighbours8() const
Definition: PFRecHit.h:90
float energy() const
rechit energy
Definition: PFRecHit.h:114
void setTime(double time)
Definition: PFRecHit.h:79
unsigned int size() const
Definition: PFRecHit.h:48
CaloCellGeometry::RepCorners RepCorners
Definition: PFRecHit.h:37
Layer
layer definition
Definition: PFLayer.h:31
bool operator>=(const PFRecHit &rhs) const
comparison >= operator
Definition: PFRecHit.h:140
RepCorners REPPointVector
Definition: PFRecHit.h:38
bool operator>(const PFRecHit &rhs) const
comparison > operator
Definition: PFRecHit.h:143
unsigned int const * Pointer
Definition: PFRecHit.h:42
PFRecHit(std::shared_ptr< const CaloCellGeometry > caloCell, unsigned int detId, PFLayer::Layer layer, float energy)
Definition: PFRecHit.h:57
CornersVec const & getCorners() const
Returns the corner points of this cell&#39;s volume.
~PFRecHit()=default
destructor
fixed size matrix
void setEnergy(float energy)
Definition: PFRecHit.h:74
CaloCellGeometry const & caloCell() const
calo cell
Definition: PFRecHit.h:104
std::array< RhoEtaPhi, k_cornerSize > RepCorners
Neighbours(Pointer ib, unsigned int n)
Definition: PFRecHit.h:44
const std::vector< unsigned short > & neighbourInfos()
Definition: PFRecHit.h:98
Neighbours neighbours() const
Definition: PFRecHit.h:94
const BasicVectorType & basicVector() const
Definition: PV3DBase.h:56
Neighbours neighbours4() const
Definition: PFRecHit.h:87
def move(src, dest)
Definition: eostools.py:510
ib
Definition: cuy.py:661
RhoEtaPhi const & positionREP() const
Definition: PFRecHit.h:131
T mag2() const
The vector magnitude squared. Equivalent to vec.dot(vec)
PFRecHit()
default constructor. Sets energy and position to zero
Definition: PFRecHit.h:55
CornersVec const & getCornersXYZ() const
rechit corners
Definition: PFRecHit.h:134