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 namespace reco {
22 
31  class PFRecHit {
32  public:
38 
39  struct Neighbours {
40  using Pointer = unsigned int const*;
42  Neighbours(Pointer ib, unsigned int n) : b(ib), e(ib + n) {}
44  Pointer begin() const { return b; }
45  Pointer end() const { return e; }
46  unsigned int size() const { return e - b; }
47  };
48 
49  enum { NONE = 0 };
51  PFRecHit() {}
52 
53  PFRecHit(std::shared_ptr<const CaloCellGeometry> caloCell,
54  unsigned int detId,
56  float energy,
57  uint32_t flags = 0)
59 
61  PFRecHit(const PFRecHit& other) = default;
62  PFRecHit(PFRecHit&& other) = default;
63  PFRecHit& operator=(const PFRecHit& other) = default;
64  PFRecHit& operator=(PFRecHit&& other) = default;
65 
67  ~PFRecHit() = default;
68 
69  void setEnergy(float energy) { energy_ = energy; }
70 
71  void addNeighbour(short x, short y, short z, unsigned int);
72  unsigned int getNeighbour(short x, short y, short z) const;
73  void setTime(double time) { time_ = time; }
74  void setDepth(int depth) { depth_ = depth; }
75  void clearNeighbours() {
76  neighbours_.clear();
77  neighbourInfos_.clear();
79  }
80 
83 
84  Neighbours neighbours() const { return buildNeighbours(neighbours_.size()); }
85 
86  const std::vector<unsigned short>& neighbourInfos() { return neighbourInfos_; }
87 
89  CaloCellGeometry const& caloCell() const { return *(caloCell_.get()); }
90  bool hasCaloCell() const { return (caloCell_ != nullptr); }
91 
93  unsigned detId() const { return detId_; }
94 
96  PFLayer::Layer layer() const { return layer_; }
97 
99  float energy() const { return energy_; }
100 
102  float time() const { return time_; }
103 
105  int depth() const { return depth_; }
106 
108  double pt2() const { return energy_ * energy_ * (position().perp2() / position().mag2()); }
109 
110  // Detector-dependent status flag
111  uint32_t flags() const { return flags_; }
112 
113  //
114  void setFlags(uint32_t flags) { flags_ = flags; }
115 
117  PositionType const& position() const { return caloCell().getPosition().basicVector(); }
118 
119  RhoEtaPhi const& positionREP() const { return caloCell().repPos(); }
120 
122  CornersVec const& getCornersXYZ() const { return caloCell().getCorners(); }
123 
124  RepCorners const& getCornersREP() const { return caloCell().getCornersREP(); }
125 
127  bool operator>=(const PFRecHit& rhs) const { return (energy_ >= rhs.energy_); }
128 
130  bool operator>(const PFRecHit& rhs) const { return (energy_ > rhs.energy_); }
131 
133  bool operator<=(const PFRecHit& rhs) const { return (energy_ <= rhs.energy_); }
134 
136  bool operator<(const PFRecHit& rhs) const { return (energy_ < rhs.energy_); }
137 
138  private:
139  Neighbours buildNeighbours(unsigned int n) const { return Neighbours(neighbours_.data(), n); }
140 
142  std::shared_ptr<const CaloCellGeometry> caloCell_ = nullptr;
143 
145  unsigned int detId_ = 0;
146 
149 
151  float energy_ = 0;
152 
154  float time_ = -1;
155 
157  int depth_ = 0;
158 
160  std::vector<unsigned int> neighbours_;
161  std::vector<unsigned short> neighbourInfos_;
162 
163  //Caching the neighbours4/8 per request of Lindsey
164  unsigned int neighbours4_ = 0;
165  unsigned int neighbours8_ = 0;
166 
167  // Detector-dependent hit status flag
168  uint32_t flags_ = 0;
169  };
170 
171 } // namespace reco
172 std::ostream& operator<<(std::ostream& out, const reco::PFRecHit& hit);
173 
174 #endif
reco::PFRecHit::time_
float time_
time
Definition: PFRecHit.h:154
reco::PFRecHit::NONE
Definition: PFRecHit.h:49
reco::PFRecHit::getCornersXYZ
CornersVec const & getCornersXYZ() const
rechit corners
Definition: PFRecHit.h:122
PFRecHitFwd.h
reco::PFRecHit::getNeighbour
unsigned int getNeighbour(short x, short y, short z) const
Definition: PFRecHit.cc:46
reco::PFRecHit::operator>=
bool operator>=(const PFRecHit &rhs) const
comparison >= operator
Definition: PFRecHit.h:127
dqmiodumpmetadata.n
n
Definition: dqmiodumpmetadata.py:28
reco::PFRecHit::buildNeighbours
Neighbours buildNeighbours(unsigned int n) const
Definition: PFRecHit.h:139
reco::PFRecHit::energy
float energy() const
rechit energy
Definition: PFRecHit.h:99
reco::PFRecHit::~PFRecHit
~PFRecHit()=default
destructor
reco::PFRecHit::neighbours8
Neighbours neighbours8() const
Definition: PFRecHit.h:82
reco::PFRecHit::operator<
bool operator<(const PFRecHit &rhs) const
comparison < operator
Definition: PFRecHit.h:136
reco::PFRecHit::flags
uint32_t flags() const
Definition: PFRecHit.h:111
reco::PFRecHit::neighbours
Neighbours neighbours() const
Definition: PFRecHit.h:84
Basic3DVector::perp2
T perp2() const
Squared magnitude of transverse component.
Definition: extBasic3DVector.h:119
CaloCellGeometry::getCorners
CornersVec const & getCorners() const
Returns the corner points of this cell's volume.
Definition: CaloCellGeometry.h:73
reco::PFRecHit::operator<=
bool operator<=(const PFRecHit &rhs) const
comparison <= operator
Definition: PFRecHit.h:133
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:46
reco::PFRecHit::Neighbours::e
Pointer e
Definition: PFRecHit.h:43
reco::PFRecHit::neighbourInfos
const std::vector< unsigned short > & neighbourInfos()
Definition: PFRecHit.h:86
reco::PFRecHit::PFRecHit
PFRecHit(std::shared_ptr< const CaloCellGeometry > caloCell, unsigned int detId, PFLayer::Layer layer, float energy, uint32_t flags=0)
Definition: PFRecHit.h:53
EZArrayFL< GlobalPoint >
Basic3DVector::mag2
T mag2() const
The vector magnitude squared. Equivalent to vec.dot(vec)
Definition: extBasic3DVector.h:113
reco::PFRecHit::Neighbours::size
unsigned int size() const
Definition: PFRecHit.h:46
operator<<
std::ostream & operator<<(std::ostream &out, const reco::PFRecHit &hit)
Definition: PFRecHit.cc:72
CaloRecHit.h
CaloCellGeometry::getCornersREP
RepCorners const & getCornersREP() const
Definition: CaloCellGeometry.h:77
reco::PFRecHit::depth
int depth() const
depth for segemntation
Definition: PFRecHit.h:105
PFLayer.h
reco::PFRecHit::REPPointVector
RepCorners REPPointVector
Definition: PFRecHit.h:36
reco::PFRecHit::PFRecHit
PFRecHit()
default constructor. Sets energy and position to zero
Definition: PFRecHit.h:51
RhoEtaPhi
Definition: PtEtaPhiMass.h:38
reco::PFRecHit::setDepth
void setDepth(int depth)
Definition: PFRecHit.h:74
reco::PFRecHit::caloCell
CaloCellGeometry const & caloCell() const
calo cell
Definition: PFRecHit.h:89
trackingPlots.other
other
Definition: trackingPlots.py:1464
PFLayer::NONE
Definition: PFLayer.h:34
reco::PFRecHit::setEnergy
void setEnergy(float energy)
Definition: PFRecHit.h:69
RefToBase.h
reco::PFRecHit::Neighbours::b
Pointer b
Definition: PFRecHit.h:43
reco::PFRecHit::Neighbours::end
Pointer end() const
Definition: PFRecHit.h:45
PFLayer::Layer
Layer
layer definition
Definition: PFLayer.h:29
CaloCellGeometry::getPosition
virtual const GlobalPoint & getPosition() const
Returns the position of reference for this cell.
Definition: CaloCellGeometry.h:80
reco::PFRecHit::addNeighbour
void addNeighbour(short x, short y, short z, unsigned int)
Definition: PFRecHit.cc:5
reco::PFRecHit::energy_
float energy_
rechit energy
Definition: PFRecHit.h:151
CaloCellGeometry
Definition: CaloCellGeometry.h:50
reco::PFRecHit::neighbours4_
unsigned int neighbours4_
Definition: PFRecHit.h:164
reco::PFRecHit::neighbourInfos_
std::vector< unsigned short > neighbourInfos_
Definition: PFRecHit.h:161
reco::PFRecHit::getCornersREP
RepCorners const & getCornersREP() const
Definition: PFRecHit.h:124
reco::PFRecHit::hasCaloCell
bool hasCaloCell() const
Definition: PFRecHit.h:90
cuy.ib
ib
Definition: cuy.py:661
reco::PFRecHit::layer_
PFLayer::Layer layer_
rechit layer
Definition: PFRecHit.h:148
PV3DBase::basicVector
const BasicVectorType & basicVector() const
Definition: PV3DBase.h:53
reco::PFRecHit::Neighbours
Definition: PFRecHit.h:39
reco::PFRecHit::time
float time() const
timing for cleaned hits
Definition: PFRecHit.h:102
reco::PFRecHit::neighbours4
Neighbours neighbours4() const
Definition: PFRecHit.h:81
reco::PFRecHit::detId
unsigned detId() const
rechit detId
Definition: PFRecHit.h:93
reco::PFRecHit::caloCell_
std::shared_ptr< const CaloCellGeometry > caloCell_
cell geometry
Definition: PFRecHit.h:142
CaloCellGeometry.h
CaloCellGeometry::CornersVec
EZArrayFL< GlobalPoint > CornersVec
Definition: CaloCellGeometry.h:57
reco::PFRecHit::depth_
int depth_
depth
Definition: PFRecHit.h:157
eostools.move
def move(src, dest)
Definition: eostools.py:511
std
Definition: JetResolutionObject.h:76
reco::PFRecHit::neighbours_
std::vector< unsigned int > neighbours_
indices to existing neighbours (1 common side)
Definition: PFRecHit.h:160
reco::PFRecHit::neighbours8_
unsigned int neighbours8_
Definition: PFRecHit.h:165
CaloCellGeometry::RepCorners
std::array< RhoEtaPhi, k_cornerSize > RepCorners
Definition: CaloCellGeometry.h:66
reco::PFRecHit::operator>
bool operator>(const PFRecHit &rhs) const
comparison > operator
Definition: PFRecHit.h:130
reco::PFRecHit::operator=
PFRecHit & operator=(const PFRecHit &other)=default
Point3D.h
CaloCellGeometry::repPos
RhoEtaPhi const & repPos() const
Definition: CaloCellGeometry.h:86
reco::PFRecHit::Neighbours::Pointer
unsigned int const * Pointer
Definition: PFRecHit.h:40
reco::PFRecHit::position
PositionType const & position() const
rechit cell centre x, y, z
Definition: PFRecHit.h:117
reco::PFRecHit::pt2
double pt2() const
rechit momentum transverse to the beam, squared.
Definition: PFRecHit.h:108
reco::PFRecHit::Neighbours::Neighbours
Neighbours(Pointer ib, unsigned int n)
Definition: PFRecHit.h:42
reco::PFRecHit::clearNeighbours
void clearNeighbours()
Definition: PFRecHit.h:75
reco::PFRecHit
Particle flow rechit (rechit + geometry and topology information). See clustering algorithm in PFClus...
Definition: PFRecHit.h:31
MillePedeFileConverter_cfg.out
out
Definition: MillePedeFileConverter_cfg.py:31
reco::PFRecHit::detId_
unsigned int detId_
cell detid
Definition: PFRecHit.h:145
reco::PFRecHit::layer
PFLayer::Layer layer() const
rechit layer
Definition: PFRecHit.h:96
reco::PFRecHit::flags_
uint32_t flags_
Definition: PFRecHit.h:168
reco::PFRecHit::Neighbours::begin
Pointer begin() const
Definition: PFRecHit.h:44
reco::PFRecHit::setTime
void setTime(double time)
Definition: PFRecHit.h:73
reco::PFRecHit::Neighbours::Neighbours
Neighbours()
Definition: PFRecHit.h:41
Vector3D.h
reco::PFRecHit::positionREP
RhoEtaPhi const & positionREP() const
Definition: PFRecHit.h:119
Basic3DVector< float >
reco::PFRecHit::setFlags
void setFlags(uint32_t flags)
Definition: PFRecHit.h:114
reco::PFRecHit::RepCorners
CaloCellGeometry::RepCorners RepCorners
Definition: PFRecHit.h:35
hit
Definition: SiStripHitEffFromCalibTree.cc:88
Point3DBase< float, GlobalTag >::BasicVectorType
BaseClass::BasicVectorType BasicVectorType
Definition: Point3DBase.h:17