CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFRecHit.cc
Go to the documentation of this file.
2 
3 using namespace reco;
4 using namespace std;
5 
6 const unsigned PFRecHit::nNeighbours_ = 8;
7 const unsigned PFRecHit::nCorners_ = 4;
8 
10  detId_(0),
12  energy_(0.),
13  rescale_(1.),
14  energyUp_(0.),
15  // seedState_(-1),
16  position_(math::XYZPoint(0.,0.,0.)),
17  posrep_(REPPoint(0.,0.,0.)) {
18 
19  cornersxyz_.reserve( nCorners_ );
20  for(unsigned i=0; i<nCorners_; i++) {
21  cornersxyz_.push_back( position_ );
22  }
23 
24  cornersrep_.reserve( nCorners_ );
25  for ( unsigned i=0; i<nCorners_; ++i ) {
26  cornersrep_.push_back(
27  REPPoint(cornersxyz_[i].Rho(),
28  cornersxyz_[i].Eta(),
29  cornersxyz_[i].Phi() ) );
30  }
31 
32 }
33 
34 
35 PFRecHit::PFRecHit(unsigned detId,
36  PFLayer::Layer layer,
37  double energy,
38  const math::XYZPoint& position,
39  const math::XYZVector& axisxyz,
40  const vector< math::XYZPoint >& cornersxyz) :
41  detId_(detId),
42  layer_(layer),
43  energy_(energy),
44  rescale_(1.),
45  energyUp_(0.),
46  // seedState_(-1),
47  position_(position),
48  posrep_( position.Rho(), position.Eta(), position.Phi() ),
49  axisxyz_(axisxyz),
50  cornersxyz_(cornersxyz)
51 {
52  cornersrep_.reserve( nCorners_ );
53  for ( unsigned i=0; i<nCorners_; ++i ) {
54  cornersrep_.push_back(
55  REPPoint(cornersxyz_[i].Rho(),
56  cornersxyz_[i].Eta(),
57  cornersxyz_[i].Phi() ) );
58  }
59 }
60 
61 PFRecHit::PFRecHit(unsigned detId,
62  PFLayer::Layer layer,
63  double energy,
64  double posx, double posy, double posz,
65  double axisx, double axisy, double axisz) :
66 
67  detId_(detId),
68  layer_(layer),
69  energy_(energy),
70  rescale_(1.),
71  energyUp_(0.),
72  // seedState_(-1),
73  position_(posx, posy, posz),
74  posrep_( position_.Rho(),
75  position_.Eta(),
76  position_.Phi() ),
77  axisxyz_(axisx, axisy, axisz) {
78 
79 
80  cornersxyz_.reserve( nCorners_ );
81  for(unsigned i=0; i<nCorners_; i++) {
82  cornersxyz_.push_back( position_ );
83  }
84 
85  cornersrep_.reserve( nCorners_ );
86  for ( unsigned i=0; i<nCorners_; ++i ) {
87  cornersrep_.push_back(
88  REPPoint(cornersxyz_[i].Rho(),
89  cornersxyz_[i].Eta(),
90  cornersxyz_[i].Phi() ) );
91  }
92 
93 }
94 
95 
97  detId_(other.detId_),
98  layer_(other.layer_),
99  energy_(other.energy_),
100  rescale_(other.rescale_),
101  energyUp_(other.energyUp_),
102  // seedState_(other.seedState_),
103  position_(other.position_),
104  posrep_(other.posrep_),
105  axisxyz_(other.axisxyz_),
106  cornersxyz_(other.cornersxyz_),
107  cornersrep_(other.cornersrep_),
108  // neighbours_(other.neighbours_),
109  neighbours4_(other.neighbours4_),
110  neighbours8_(other.neighbours8_),
111  neighboursIds4_(other.neighboursIds4_),
112  neighboursIds8_(other.neighboursIds8_)
113 {}
114 
115 
117 {}
118 
119 
120 const PFRecHit::REPPoint&
122  // if( posrep_ == REPPoint() ) {
123  // posrep_.SetCoordinates( position_.Rho(), position_.Eta(), position_.Phi() );
124  // }
125  return posrep_;
126 }
127 
128 
129 void
131  posrep_.SetCoordinates( position_.Rho(), position_.Eta(), position_.Phi() );
132  for ( unsigned i=0; i<nCorners_; ++i ) {
133  cornersrep_[i].SetCoordinates(cornersxyz_[i].Rho(),cornersxyz_[i].Eta(),cornersxyz_[i].Phi() );
134  }
135 }
136 
137 
138 // void PFRecHit::setNeighbours( const vector< unsigned >& neighbours ) {
139 // if( neighbours.size() != nNeighbours_ )
140 // throw cms::Exception("CellNeighbourVector")
141 // << "number of neighbours must be nNeighbours_";
142 
143 // // neighbours_.clear();
144 // // neighboursIds4_.clear();
145 // // neighboursIds8_.clear();
146 // neighbours4_.clear();
147 // neighbours8_.clear();
148 
149 // // neighboursIds4_.reserve( PFRecHit::nNeighbours_ );
150 // // neighboursIds8_.reserve( PFRecHit::nNeighbours_ );
151 
152 // // space is reserved, but this vectors will not always have size 4 or 8.
153 // // they contain the indices to the neighbours that are present in the rechit
154 // // collection
155 // neighbours4_.reserve( PFRecHit::nNeighbours_ );
156 // neighbours8_.reserve( PFRecHit::nNeighbours_ );
157 
158 // // neighboursIds_.reserve(nNeighbours_);
159 // // neighbours_.reserve(nNeighbours_);
160 
161 // // neighbours_ = neighbours;
162 
163 // for(unsigned i=0; i<neighbours.size(); i++) {
164 // if( neighbours[i] ) {
165 // neighbours8_.push_back( neighbours[i] );
166 // // neighboursIds8_.push_back( neighbours[i]->detId() );
167 // if( !(i%2) ) {
168 // neighbours4_.push_back( neighbours[i] );
169 // // neighboursIds4_.push_back( neighbours[i]->detId() );
170 // }
171 // }
172 // }
173 // }
174 
175 
177  neighbours4_.push_back( index );
178  neighbours8_.push_back( index );
179 }
180 
182  neighbours8_.push_back( index );
183 }
184 
185 
186 
187 
188 // void PFRecHit::findPtrsToNeighbours( const std::map<unsigned,
189 // reco::PFRecHit* >& allhits ) {
190 
191 // neighbours4_.clear();
192 // neighbours8_.clear();
193 
194 // typedef std::map<unsigned, reco::PFRecHit* >::const_iterator IDH;
195 
196 // for(unsigned inid = 0; inid<neighboursIds8_.size(); inid++) {
197 
198 // IDH ineighbour = allhits.find( neighboursIds8_[inid] );
199 // if( ineighbour != allhits.end() ) {
200 // neighbours8_.push_back( ineighbour->second );
201 // }
202 // }
203 
204 
205 // for(unsigned inid = 0; inid<neighboursIds4_.size(); inid++) {
206 
207 // IDH ineighbour = allhits.find( neighboursIds4_[inid] );
208 // if( ineighbour != allhits.end() ) {
209 // neighbours4_.push_back( ineighbour->second );
210 // }
211 // }
212 
213 // }
214 
215 
216 void PFRecHit::setNWCorner( double posx, double posy, double posz ) {
217  setCorner(0, posx, posy, posz);
218 }
219 
220 
221 void PFRecHit::setSWCorner( double posx, double posy, double posz ) {
222  setCorner(1, posx, posy, posz);
223 }
224 
225 
226 void PFRecHit::setSECorner( double posx, double posy, double posz ) {
227  setCorner(2, posx, posy, posz);
228 }
229 
230 
231 void PFRecHit::setNECorner( double posx, double posy, double posz ) {
232  setCorner(3, posx, posy, posz);
233 }
234 
235 
236 void PFRecHit::setCorner( unsigned i, double posx, double posy, double posz ) {
237  assert( cornersxyz_.size() == nCorners_);
238  assert( cornersrep_.size() == nCorners_);
239  assert( i<cornersxyz_.size() );
240 
241  cornersxyz_[i] = math::XYZPoint( posx, posy, posz);
242  cornersrep_[i] = REPPoint(cornersxyz_[i].Rho(),
243  cornersxyz_[i].Eta(),
244  cornersxyz_[i].Phi() );
245 }
246 
247 
248 bool PFRecHit::isNeighbour4(unsigned id) const {
249 
250  for(unsigned i=0; i<neighbours4_.size(); i++ )
251  if( id == neighbours4_[i] ) return true;
252 
253  return false;
254 }
255 
256 
257 bool PFRecHit::isNeighbour8(unsigned id) const {
258 
259  for(unsigned i=0; i<neighbours8_.size(); i++ )
260  if( id == neighbours8_[i] ) return true;
261 
262  return false;
263 }
264 
265 
266 void PFRecHit::size(double& deta, double& dphi) const {
267 
268  double minphi=9999;
269  double maxphi=-9999;
270  double mineta=9999;
271  double maxeta=-9999;
272  for ( unsigned ic=0; ic<cornersrep_.size(); ++ic ) {
273  double eta = cornersrep_[ic].Eta();
274  double phi = cornersrep_[ic].Phi();
275 
276  if(phi>maxphi) maxphi=phi;
277  if(phi<minphi) minphi=phi;
278  if(eta>maxeta) maxeta=eta;
279  if(eta<mineta) mineta=eta;
280  }
281 
282  deta = maxeta - mineta;
283  dphi = maxphi - minphi;
284 }
285 
286 
287 ostream& reco::operator<<(ostream& out, const reco::PFRecHit& hit) {
288 
289  if(!out) return out;
290 
291  // reco::PFRecHit& nshit = const_cast<reco::PFRecHit& >(hit);
292  // const reco::PFRecHit::REPPoint& posrep = nshit.positionREP();
293 
294  const math::XYZPoint& posxyz = hit.position();
295 
296  out<<"hit id:"<<hit.detId()
297  <<" l:"<<hit.layer()
298  <<" E:"<<hit.energy()
299  <<" t:"<<hit.time()
300  <<" rep:"<<posxyz.Rho()<<","<<posxyz.Eta()<<","<<posxyz.Phi()<<"| N:";
301  // <<" SEED: "<<hit.seedState_;
302  for(unsigned i=0; i<hit.neighbours8_.size(); i++ ) {
303  out<<hit.neighbours8_[i]<<" ";
304  }
305 
306 
307  // out<<endl;
308  // out<<"neighbours "<<endl;
309  // for(unsigned i=0; i<hit.neighbours8_.size(); i++ ) {
310  // out<<"\t"<< hit.neighbours8_[i]->detId()<<endl;
311  // }
312  // out<<"--"<<endl;
313  // for(unsigned i=0; i<hit.neighboursIds8_.size(); i++ ) {
314  // out<<"\t"<< (hit.neighboursIds8_[i])<<endl;
315  // }
316 
317  // bool printcorners = false;
318 
319  // if(printcorners) {
320 
321  // out<<endl<<"corners : "<<endl;
322  // const std::vector< math::XYZPoint >& corners = hit.getCornersXYZ();
323  // for(unsigned i=0; i<corners.size(); i++) {
324  // out<<"\t"<<corners[i].X()<<","<<corners[i].Y()<<","<<corners[i].Z()<<endl;
325  // }
326  // }
327 
328  return out;
329 }
void setSECorner(double posx, double posy, double posz)
Definition: PFRecHit.cc:226
void add4Neighbour(unsigned index)
Definition: PFRecHit.cc:176
int i
Definition: DBlmapReader.cc:9
static const char layer_[]
std::vector< unsigned > neighbours4_
indices to existing neighbours (1 common side)
Definition: PFRecHit.h:235
REPPoint posrep_
rechit cell centre: rho, eta, phi (transient)
Definition: PFRecHit.h:225
ROOT::Math::PositionVector3D< ROOT::Math::CylindricalEta3D< Double32_t > > REPPoint
Definition: PFRecHit.h:39
std::vector< math::XYZPoint > cornersxyz_
rechit cell corners
Definition: PFRecHit.h:231
void add8Neighbour(unsigned index)
Definition: PFRecHit.cc:181
void setNECorner(double posx, double posy, double posz)
Definition: PFRecHit.cc:231
unsigned detId() const
rechit detId
Definition: PFRecHit.h:105
T eta() const
static const unsigned nNeighbours_
number of neighbours
Definition: PFRecHit.h:247
static int position[TOTALCHAMBERS][3]
Definition: ReadPGInfo.cc:509
PFLayer::Layer layer() const
rechit layer
Definition: PFRecHit.h:108
virtual ~PFRecHit()
destructor
Definition: PFRecHit.cc:116
void setSWCorner(double posx, double posy, double posz)
Definition: PFRecHit.cc:221
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:31
layer definition for PFRecHit and PFCluster
Definition: PFLayer.h:21
void size(double &deta, double &dphi) const
Definition: PFRecHit.cc:266
std::vector< unsigned > neighbours8_
indices to existing neighbours (1 common side or diagonal)
Definition: PFRecHit.h:238
void setNWCorner(double posx, double posy, double posz)
search for pointers to neighbours, using neighbours&#39; DetId.
Definition: PFRecHit.cc:216
REPPointVector cornersrep_
Definition: PFRecHit.h:232
tuple out
Definition: dbtoconf.py:99
Layer
layer definition
Definition: PFLayer.h:31
const math::XYZPoint & position() const
is seed ?
Definition: PFRecHit.h:141
void calculatePositionREP()
calculates rho eta phi position once and for all
Definition: PFRecHit.cc:130
math::XYZPoint position_
is this a seed ? (-1:unknown, 0:no, 1 yes) (transient)
Definition: PFRecHit.h:222
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
bool isNeighbour4(unsigned id) const
Definition: PFRecHit.cc:248
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
double energy() const
rechit energy
Definition: PFRecHit.h:111
bool isNeighbour8(unsigned id) const
Definition: PFRecHit.cc:257
static const unsigned nCorners_
number of corners
Definition: PFRecHit.h:250
double time() const
timing for cleaned hits
Definition: PFRecHit.h:117
void setCorner(unsigned i, double posx, double posy, double posz)
set position of one of the corners
Definition: PFRecHit.cc:236
PFRecHit()
default constructor. Sets energy and position to zero
Definition: PFRecHit.cc:9
const REPPoint & positionREP() const
rechit cell centre rho, eta, phi. call calculatePositionREP before !
Definition: PFRecHit.cc:121
Definition: DDAxes.h:10