CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 <iostream>
8 
10 #include "Rtypes.h"
12 #include "Math/GenVector/PositionVector3D.h"
15 
16 //C decide what is the default rechit index.
17 //C maybe 0 ? -> compression
18 //C then the position is index-1.
19 //C provide a helper class to access the rechit.
20 
23 
24 
25 namespace reco {
26 
35  class PFRecHit {
36 
37  public:
38 
39  // Next typedef uses double in ROOT 6 rather than Double32_t due to a bug in ROOT 5,
40  // which otherwise would make ROOT5 files unreadable in ROOT6. This does not increase
41  // the size on disk, because due to the bug, double was actually stored on disk in ROOT 5.
42  typedef ROOT::Math::PositionVector3D<ROOT::Math::CylindricalEta3D<double> > REPPoint;
43 
44  typedef std::vector<REPPoint> REPPointVector;
45 
46  enum {
47  NONE=0
48  };
50  PFRecHit();
51 
53  PFRecHit(unsigned detId,
55  double energy,
56  const math::XYZPoint& posxyz,
57  const math::XYZVector& axisxyz,
58  const std::vector< math::XYZPoint >& cornersxyz);
59 
60  PFRecHit(unsigned detId,
62  double energy,
63  double posx, double posy, double posz,
64  double axisx, double axisy, double axisz);
65 
67  PFRecHit(const PFRecHit& other);
68 
70  virtual ~PFRecHit();
71 
72  void setEnergy( double energy) { energy_ = energy; }
73 
74  void calculatePositionREP();
75 
76  void addNeighbour(short x,short y, short z,const PFRecHitRef&);
77  const PFRecHitRef getNeighbour(short x,short y, short z);
78  void setTime( double time) { time_ = time; }
79  void setDepth( int depth) { depth_ = depth; }
80  void clearNeighbours() {
82  }
83 
84  const PFRecHitRefVector& neighbours4() const {
85  return neighbours4_;
86  }
87  const PFRecHitRefVector& neighbours8() const {
88  return neighbours8_;
89  }
90 
91  const PFRecHitRefVector& neighbours() const {
92  return neighbours_;
93  }
94 
95  const std::vector<unsigned short>& neighbourInfos() {
96  return neighbourInfos_;
97  }
98 
99 
100  void setNWCorner( double posx, double posy, double posz );
101  void setSWCorner( double posx, double posy, double posz );
102  void setSECorner( double posx, double posy, double posz );
103  void setNECorner( double posx, double posy, double posz );
104 
106  unsigned detId() const {return detId_;}
107 
109  PFLayer::Layer layer() const { return layer_; }
110 
112  double energy() const { return energy_; }
113 
114 
116  double time() const { return time_; }
117 
119  int depth() const { return depth_; }
120 
122  double pt2() const { return energy_ * energy_ *
123  ( position_.X()*position_.X() +
124  position_.Y()*position_.Y() ) /
125  ( position_.X()*position_.X() +
126  position_.Y()*position_.Y() +
127  position_.Z()*position_.Z()) ; }
128 
129 
131  const math::XYZPoint& position() const { return position_; }
132 
133  const REPPoint& positionREP() const { return positionrep_; }
134 
135 
137  const math::XYZVector& getAxisXYZ() const { return axisxyz_; }
138 
140  const std::vector< math::XYZPoint >& getCornersXYZ() const
141  { return cornersxyz_; }
142 
143  const std::vector<REPPoint>& getCornersREP() const { return cornersrep_; }
144 
145  void size(double& deta, double& dphi) const;
146 
148  bool operator>=(const PFRecHit& rhs) const { return (energy_>=rhs.energy_); }
149 
151  bool operator> (const PFRecHit& rhs) const { return (energy_> rhs.energy_); }
152 
154  bool operator<=(const PFRecHit& rhs) const { return (energy_<=rhs.energy_); }
155 
157  bool operator< (const PFRecHit& rhs) const { return (energy_< rhs.energy_); }
158 
159  friend std::ostream& operator<<(std::ostream& out,
160  const reco::PFRecHit& hit);
161 
163  return originalRecHit_;
164  }
165 
166  template<typename T>
167  void setOriginalRecHit(const T& rh) {
169  }
170 
171  private:
172 
173  // original rechit
175 
177  unsigned detId_;
178 
181 
183  double energy_;
184 
186  double time_;
187 
188 
190  int depth_;
191 
194 
197 
200 
202  std::vector< math::XYZPoint > cornersxyz_;
203 
205  std::vector< REPPoint > cornersrep_;
206 
209  std::vector< unsigned short > neighbourInfos_;
210 
211  //Caching the neighbours4/8 per request of Lindsey
214 
215 
217  static const unsigned nCorners_;
218 
220  void setCorner( unsigned i, double posx, double posy, double posz );
221  };
222 
223 }
224 
225 #endif
void setSECorner(double posx, double posy, double posz)
Definition: PFRecHit.cc:100
void setDepth(int depth)
Definition: PFRecHit.h:79
int i
Definition: DBlmapReader.cc:9
std::vector< REPPoint > REPPointVector
Definition: PFRecHit.h:44
std::vector< math::XYZPoint > cornersxyz_
rechit cell corners
Definition: PFRecHit.h:202
friend std::ostream & operator<<(std::ostream &out, const reco::PFRecHit &hit)
void setNECorner(double posx, double posy, double posz)
Definition: PFRecHit.cc:105
double pt2() const
rechit momentum transverse to the beam, squared.
Definition: PFRecHit.h:122
const PFRecHitRefVector & neighbours4() const
Definition: PFRecHit.h:84
unsigned detId_
C cell detid - should be detid or index in collection ?
Definition: PFRecHit.h:177
PFRecHitRefVector neighbours4_
Definition: PFRecHit.h:212
unsigned detId() const
rechit detId
Definition: PFRecHit.h:106
int depth_
depth
Definition: PFRecHit.h:190
double time_
time
Definition: PFRecHit.h:186
void addNeighbour(short x, short y, short z, const PFRecHitRef &)
Definition: PFRecHit.cc:133
PFLayer::Layer layer_
rechit layer
Definition: PFRecHit.h:180
std::vector< unsigned short > neighbourInfos_
Definition: PFRecHit.h:209
bool operator<=(const PFRecHit &rhs) const
comparison &lt;= operator
Definition: PFRecHit.h:154
bool operator<(const PFRecHit &rhs) const
comparison &lt; operator
Definition: PFRecHit.h:157
double energy_
rechit energy
Definition: PFRecHit.h:183
int depth() const
depth for segemntation
Definition: PFRecHit.h:119
PFLayer::Layer layer() const
rechit layer
Definition: PFRecHit.h:109
virtual ~PFRecHit()
destructor
Definition: PFRecHit.cc:86
void setSWCorner(double posx, double posy, double posz)
Definition: PFRecHit.cc:95
const std::vector< math::XYZPoint > & getCornersXYZ() const
rechit corners
Definition: PFRecHit.h:140
void setEnergy(double energy)
Definition: PFRecHit.h:72
const PFRecHitRef getNeighbour(short x, short y, short z)
Definition: PFRecHit.cc:171
void clearNeighbours()
Definition: PFRecHit.h:80
Particle flow rechit (rechit + geometry and topology information). See clustering algorithm in PFClus...
Definition: PFRecHit.h:35
math::XYZVector axisxyz_
rechit cell axisxyz
Definition: PFRecHit.h:199
const math::XYZVector & getAxisXYZ() const
rechit cell axis x, y, z
Definition: PFRecHit.h:137
const std::vector< REPPoint > & getCornersREP() const
Definition: PFRecHit.h:143
const PFRecHitRefVector & neighbours8() const
Definition: PFRecHit.h:87
void size(double &deta, double &dphi) const
Definition: PFRecHit.cc:196
edm::RefToBase< CaloRecHit > originalRecHit_
Definition: PFRecHit.h:174
void setTime(double time)
Definition: PFRecHit.h:78
void setNWCorner(double posx, double posy, double posz)
Definition: PFRecHit.cc:90
const edm::RefToBase< CaloRecHit > & originalRecHit() const
Definition: PFRecHit.h:162
tuple out
Definition: dbtoconf.py:99
Layer
layer definition
Definition: PFLayer.h:31
const math::XYZPoint & position() const
rechit cell centre x, y, z
Definition: PFRecHit.h:131
void calculatePositionREP()
Definition: PFRecHit.cc:121
ROOT::Math::PositionVector3D< ROOT::Math::CylindricalEta3D< double > > REPPoint
Definition: PFRecHit.h:42
PFRecHitRefVector neighbours8_
Definition: PFRecHit.h:213
bool operator>=(const PFRecHit &rhs) const
comparison &gt;= operator
Definition: PFRecHit.h:148
math::XYZPoint position_
rechit cell centre: x, y, z
Definition: PFRecHit.h:193
void clear()
Clear the vector.
Definition: RefVector.h:139
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
bool operator>(const PFRecHit &rhs) const
comparison &gt; operator
Definition: PFRecHit.h:151
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
PFRecHitRefVector neighbours_
indices to existing neighbours (1 common side)
Definition: PFRecHit.h:208
std::vector< REPPoint > cornersrep_
rechit cell corners rho/eta/phi
Definition: PFRecHit.h:205
double energy() const
rechit energy
Definition: PFRecHit.h:112
REPPoint positionrep_
rechit cell centre: rho, eta, phi (transient)
Definition: PFRecHit.h:196
const REPPoint & positionREP() const
Definition: PFRecHit.h:133
static const unsigned nCorners_
number of corners
Definition: PFRecHit.h:217
double time() const
timing for cleaned hits
Definition: PFRecHit.h:116
const std::vector< unsigned short > & neighbourInfos()
Definition: PFRecHit.h:95
long double T
void setCorner(unsigned i, double posx, double posy, double posz)
set position of one of the corners
Definition: PFRecHit.cc:110
const PFRecHitRefVector & neighbours() const
Definition: PFRecHit.h:91
PFRecHit()
default constructor. Sets energy and position to zero
Definition: PFRecHit.cc:7
void setOriginalRecHit(const T &rh)
Definition: PFRecHit.h:167