CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TrackUtils.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: Tracks
4 // Class : TrackUtils
5 //
6 
7 // system include files
8 #include "TEveTrack.h"
9 #include "TEveTrackPropagator.h"
10 #include "TEveStraightLineSet.h"
11 #include "TEveVSDStructs.h"
12 #include "TEveGeoNode.h"
13 
14 // user include files
20 
22 
24 
26 
42 
47 
51 
53 
56 
58 
59 namespace fireworks {
60 
61 static const double MICRON = 1./1000./10.;
62 static const double PITCHX = 100*MICRON;
63 static const double PITCHY = 150*MICRON;
64 static const int BIG_PIX_PER_ROC_Y = 2; // in y direction, cols
65 static const int COLS_PER_ROC = 52; // Num of Rows per ROC
66 static const int ROWS_PER_ROC = 80; // Num of cols per ROC
67 static const int BIG_PIX_PER_ROC_X = 1; // in x direction, rows
68 
69 // -- Si module names for printout
70 static const std::string subdets[7] = { "UNKNOWN", "PXB", "PXF", "TIB", "TID", "TOB", "TEC" };
71 
72 TEveTrack*
73 prepareTrack( const reco::Track& track,
74  TEveTrackPropagator* propagator,
75  const std::vector<TEveVector>& extraRefPoints )
76 {
77  // To make use of all available information, we have to order states
78  // properly first. Propagator should take care of y=0 transition.
79 
80  std::vector<State> refStates;
81  TEveVector trackMomentum(track.px(), track.py(), track.pz());
82  refStates.push_back(State(TEveVector(track.vx(), track.vy(), track.vz()),
83  trackMomentum));
84  if( track.extra().isAvailable() ) {
85  if (track.innerOk()) {
86  const reco::TrackBase::Point &v = track.innerPosition();
87  const reco::TrackBase::Vector &p = track.innerMomentum();
88  refStates.push_back(State(TEveVector(v.x(), v.y(), v.z()), TEveVector(p.x(), p.y(), p.z())));
89  }
90  if (track.outerOk()) {
91  const reco::TrackBase::Point &v = track.outerPosition();
92  const reco::TrackBase::Vector &p = track.outerMomentum();
93  refStates.push_back(State(TEveVector(v.x(), v.y(), v.z()), TEveVector(p.x(), p.y(), p.z())));
94  }
95  }
96  for( std::vector<TEveVector>::const_iterator point = extraRefPoints.begin(), pointEnd = extraRefPoints.end();
97  point != pointEnd; ++point )
98  refStates.push_back(State(*point));
99  if( track.pt()>1 )
100  std::sort( refStates.begin(), refStates.end(), StateOrdering(trackMomentum) );
101 
102  // * if the first state has non-zero momentum use it as a starting point
103  // and all other points as PathMarks to follow
104  // * if the first state has only position, try the last state. If it has
105  // momentum we propagate backword, if not, we look for the first one
106  // on left that has momentum and ignore all earlier.
107  //
108 
109  TEveRecTrack t;
110  t.fBeta = 1.;
111  t.fSign = track.charge();
112 
113  if( refStates.front().valid ) {
114  t.fV = refStates.front().position;
115  t.fP = refStates.front().momentum;
116  TEveTrack* trk = new TEveTrack( &t, propagator );
117  for( unsigned int i(1); i<refStates.size()-1; ++i) {
118  if( refStates[i].valid )
119  trk->AddPathMark( TEvePathMark( TEvePathMark::kReference, refStates[i].position, refStates[i].momentum ) );
120  else
121  trk->AddPathMark( TEvePathMark( TEvePathMark::kDaughter, refStates[i].position ) );
122  }
123  if( refStates.size()>1 ) {
124  trk->AddPathMark( TEvePathMark( TEvePathMark::kDecay, refStates.back().position ) );
125  }
126  return trk;
127  }
128 
129  if( refStates.back().valid ) {
130  t.fSign = (-1)*track.charge();
131  t.fV = refStates.back().position;
132  t.fP = refStates.back().momentum * (-1.0f);
133  TEveTrack* trk = new TEveTrack( &t, propagator );
134  unsigned int i( refStates.size()-1 );
135  for(; i>0; --i) {
136  if ( refStates[i].valid )
137  trk->AddPathMark( TEvePathMark( TEvePathMark::kReference, refStates[i].position, refStates[i].momentum*(-1.0f) ) );
138  else
139  trk->AddPathMark( TEvePathMark( TEvePathMark::kDaughter, refStates[i].position ) );
140  }
141  if ( refStates.size()>1 ) {
142  trk->AddPathMark( TEvePathMark( TEvePathMark::kDecay, refStates.front().position ) );
143  }
144  return trk;
145  }
146 
147  unsigned int i(0);
148  while( i<refStates.size() && !refStates[i].valid ) ++i;
149  assert( i < refStates.size() );
150 
151  t.fV = refStates[i].position;
152  t.fP = refStates[i].momentum;
153  TEveTrack* trk = new TEveTrack( &t, propagator );
154  for( unsigned int j(i+1); j<refStates.size()-1; ++j ) {
155  if( refStates[i].valid )
156  trk->AddPathMark( TEvePathMark( TEvePathMark::kReference, refStates[i].position, refStates[i].momentum ) );
157  else
158  trk->AddPathMark( TEvePathMark( TEvePathMark::kDaughter, refStates[i].position ) );
159  }
160  if ( i < refStates.size() ) {
161  trk->AddPathMark( TEvePathMark( TEvePathMark::kDecay, refStates.back().position ) );
162  }
163  return trk;
164 }
165 
166 // Transform measurement to local coordinates in X dimension
167 float
168 pixelLocalX( const double mpx, const int nrows )
169 {
170  // Calculate the edge of the active sensor with respect to the center,
171  // that is simply the half-size.
172  // Take into account large pixels
173  const double xoffset = -( nrows + BIG_PIX_PER_ROC_X * nrows / ROWS_PER_ROC ) / 2. * PITCHX;
174 
175  // Measurement to local transformation for X coordinate
176  // X coordinate is in the ROC row number direction
177  // Copy from RectangularPixelTopology::localX implementation
178  int binoffx = int(mpx); // truncate to int
179  double fractionX = mpx - binoffx; // find the fraction
180  double local_PITCHX = PITCHX; // defaultpitch
181  if( binoffx > 80 ) { // ROC 1 - handles x on edge cluster
182  binoffx = binoffx + 2;
183  } else if( binoffx == 80 ) { // ROC 1
184  binoffx = binoffx+1;
185  local_PITCHX = 2 * PITCHX;
186  } else if( binoffx == 79 ) { // ROC 0
187  binoffx = binoffx + 0;
188  local_PITCHX = 2 * PITCHX;
189  } else if( binoffx >= 0 ) { // ROC 0
190  binoffx = binoffx + 0;
191  }
192 
193  // The final position in local coordinates
194  double lpX = double( binoffx * PITCHX ) + fractionX * local_PITCHX + xoffset;
195 
196  return lpX;
197 }
198 
199 // Transform measurement to local coordinates in Y dimension
200 float
201 pixelLocalY( const double mpy, const int ncols )
202 {
203  // Calculate the edge of the active sensor with respect to the center,
204  // that is simply the half-size.
205  // Take into account large pixels
206  double yoffset = -( ncols + BIG_PIX_PER_ROC_Y * ncols / COLS_PER_ROC ) / 2. * PITCHY;
207 
208  // Measurement to local transformation for Y coordinate
209  // Y is in the ROC column number direction
210  // Copy from RectangularPixelTopology::localY implementation
211  int binoffy = int( mpy ); // truncate to int
212  double fractionY = mpy - binoffy; // find the fraction
213  double local_PITCHY = PITCHY; // defaultpitch
214 
215  if( binoffy>416 ) { // ROC 8, not real ROC
216  binoffy = binoffy+17;
217  } else if( binoffy == 416 ) { // ROC 8
218  binoffy = binoffy + 16;
219  local_PITCHY = 2 * PITCHY;
220  } else if( binoffy == 415 ) { // ROC 7, last big pixel
221  binoffy = binoffy + 15;
222  local_PITCHY = 2 * PITCHY;
223  } else if( binoffy > 364 ) { // ROC 7
224  binoffy = binoffy + 15;
225  } else if( binoffy == 364 ) { // ROC 7
226  binoffy = binoffy + 14;
227  local_PITCHY = 2 * PITCHY;
228  } else if( binoffy == 363 ) { // ROC 6
229  binoffy = binoffy + 13;
230  local_PITCHY = 2 * PITCHY;
231  } else if( binoffy > 312 ) { // ROC 6
232  binoffy = binoffy + 13;
233  } else if( binoffy == 312 ) { // ROC 6
234  binoffy = binoffy + 12;
235  local_PITCHY = 2 * PITCHY;
236  } else if( binoffy == 311 ) { // ROC 5
237  binoffy = binoffy + 11;
238  local_PITCHY = 2 * PITCHY;
239  } else if( binoffy > 260 ) { // ROC 5
240  binoffy = binoffy + 11;
241  } else if( binoffy == 260 ) { // ROC 5
242  binoffy = binoffy + 10;
243  local_PITCHY = 2 * PITCHY;
244  } else if( binoffy == 259 ) { // ROC 4
245  binoffy = binoffy + 9;
246  local_PITCHY = 2 * PITCHY;
247  } else if( binoffy > 208 ) { // ROC 4
248  binoffy = binoffy + 9;
249  } else if(binoffy == 208 ) { // ROC 4
250  binoffy = binoffy + 8;
251  local_PITCHY = 2 * PITCHY;
252  } else if( binoffy == 207 ) { // ROC 3
253  binoffy = binoffy + 7;
254  local_PITCHY = 2 * PITCHY;
255  } else if( binoffy > 156 ) { // ROC 3
256  binoffy = binoffy + 7;
257  } else if( binoffy == 156 ) { // ROC 3
258  binoffy = binoffy + 6;
259  local_PITCHY = 2 * PITCHY;
260  } else if( binoffy == 155 ) { // ROC 2
261  binoffy = binoffy + 5;
262  local_PITCHY = 2 * PITCHY;
263  } else if( binoffy > 104 ) { // ROC 2
264  binoffy = binoffy + 5;
265  } else if( binoffy == 104 ) { // ROC 2
266  binoffy = binoffy + 4;
267  local_PITCHY = 2 * PITCHY;
268  } else if( binoffy == 103 ) { // ROC 1
269  binoffy = binoffy + 3;
270  local_PITCHY = 2 * PITCHY;
271  } else if( binoffy > 52 ) { // ROC 1
272  binoffy = binoffy + 3;
273  } else if( binoffy == 52 ) { // ROC 1
274  binoffy = binoffy + 2;
275  local_PITCHY = 2 * PITCHY;
276  } else if( binoffy == 51 ) { // ROC 0
277  binoffy = binoffy + 1;
278  local_PITCHY = 2 * PITCHY;
279  } else if( binoffy > 0 ) { // ROC 0
280  binoffy=binoffy + 1;
281  } else if( binoffy == 0 ) { // ROC 0
282  binoffy = binoffy + 0;
283  local_PITCHY = 2 * PITCHY;
284  }
285 
286  // The final position in local coordinates
287  double lpY = double( binoffy * PITCHY ) + fractionY * local_PITCHY + yoffset;
288 
289  return lpY;
290 }
291 
292 //
293 // Returns strip geometry in local coordinates of a detunit.
294 // The strip is a line from a localTop to a localBottom point.
295 void localSiStrip( short strip, float* localTop, float* localBottom, const float* pars, unsigned int id )
296 {
297  Float_t topology = pars[0];
298  Float_t halfStripLength = pars[2] * 0.5;
299 
300  Double_t localCenter[3] = { 0.0, 0.0, 0.0 };
301  localTop[1] = halfStripLength;
302  localBottom[1] = -halfStripLength;
303 
304  if( topology == 1 ) // RadialStripTopology
305  {
306  // stripAngle = phiOfOneEdge + strip * angularWidth
307  // localY = originToIntersection * tan( stripAngle )
308  Float_t stripAngle = tan( pars[5] + strip * pars[6] );
309  Float_t delta = halfStripLength * stripAngle;
310  localCenter[0] = pars[4] * stripAngle;
311  localTop[0] = localCenter[0] + delta;
312  localBottom[0] = localCenter[0] - delta;
313  }
314  else if( topology == 2 ) // RectangularStripTopology
315  {
316  // offset = -numberOfStrips/2. * pitch
317  // localY = strip * pitch + offset
318  Float_t offset = -pars[1] * 0.5 * pars[3];
319  localCenter[0] = strip * pars[3] + offset;
320  localTop[0] = localCenter[0];
321  localBottom[0] = localCenter[0];
322  }
323  else if( topology == 3 ) // TrapezoidalStripTopology
324  {
326  << "did not expect TrapezoidalStripTopology of "
327  << id << std::endl;
328  }
329  else if( pars[0] == 0 ) // StripTopology
330  {
332  << "did not find StripTopology of "
333  << id << std::endl;
334  }
335 }
336 
337 //______________________________________________________________________________
338 
339 void
340 setupAddElement(TEveElement* el, TEveElement* parent, const FWEventItem* item, bool master, bool color)
341 {
342  if (master)
343  {
344  el->CSCTakeAnyParentAsMaster();
345  el->SetPickable(true);
346  }
347 
348  if (color)
349  {
350  el->CSCApplyMainColorToMatchingChildren();
351  el->CSCApplyMainTransparencyToMatchingChildren();
352  el->SetMainColor(item->defaultDisplayProperties().color());
354  && (item->defaultDisplayProperties().transparency() <= 100));
355  el->SetMainTransparency(item->defaultDisplayProperties().transparency());
356  }
357  parent->AddElement(el);
358 }
359 
360 //______________________________________________________________________________
361 
363 {
364  const SiStripCluster* cluster = 0;
365 
366  if( const SiStripRecHit2D* hit2D = dynamic_cast<const SiStripRecHit2D*>( rechit ))
367  {
368  fwLog( fwlog::kDebug ) << "hit 2D ";
369 
370  cluster = hit2D->cluster().get();
371  }
372  if( cluster == 0 )
373  {
374  if( const SiStripRecHit1D* hit1D = dynamic_cast<const SiStripRecHit1D*>( rechit ))
375  {
376  fwLog( fwlog::kDebug ) << "hit 1D ";
377 
378  cluster = hit1D->cluster().get();
379  }
380  }
381  return cluster;
382 }
383 
384 void
385 addSiStripClusters( const FWEventItem* iItem, const reco::Track &t, class TEveElement *tList, bool addNearbyClusters, bool master )
386 {
387  // master is true if the product is for proxy builder
388  const FWGeometry *geom = iItem->getGeom();
389 
390  const edmNew::DetSetVector<SiStripCluster> * allClusters = 0;
391  if( addNearbyClusters )
392  {
393  for( trackingRecHit_iterator it = t.recHitsBegin(), itEnd = t.recHitsEnd(); it != itEnd; ++it )
394  {
395  if( typeid( **it ) == typeid( SiStripRecHit2D ))
396  {
397  const SiStripRecHit2D &hit = static_cast<const SiStripRecHit2D &>( **it );
398  if( hit.cluster().isNonnull() && hit.cluster().isAvailable())
399  {
401  iItem->getEvent()->get(hit.cluster().id(), allClustersHandle);
402  allClusters = allClustersHandle.product();
403  break;
404  }
405  }
406  else if( typeid( **it ) == typeid( SiStripRecHit1D ))
407  {
408  const SiStripRecHit1D &hit = static_cast<const SiStripRecHit1D &>( **it );
409  if( hit.cluster().isNonnull() && hit.cluster().isAvailable())
410  {
412  iItem->getEvent()->get(hit.cluster().id(), allClustersHandle);
413  allClusters = allClustersHandle.product();
414  break;
415  }
416  }
417  }
418  }
419 
420  for( trackingRecHit_iterator it = t.recHitsBegin(), itEnd = t.recHitsEnd(); it != itEnd; ++it )
421  {
422  unsigned int rawid = (*it)->geographicalId();
423  if( ! geom->contains( rawid ))
424  {
426  << "failed to get geometry of SiStripCluster with detid: "
427  << rawid << std::endl;
428 
429  continue;
430  }
431 
432  const float* pars = geom->getParameters( rawid );
433 
434  // -- get phi from SiStripHit
435  auto rechitRef = *it;
436  const TrackingRecHit *rechit = &( *rechitRef );
437  const SiStripCluster *cluster = extractClusterFromTrackingRecHit( rechit );
438 
439  if( cluster )
440  {
441  if( allClusters != 0 )
442  {
443  const edmNew::DetSet<SiStripCluster> & clustersOnThisDet = (*allClusters)[rechit->geographicalId().rawId()];
444 
445  for( edmNew::DetSet<SiStripCluster>::const_iterator itc = clustersOnThisDet.begin(), edc = clustersOnThisDet.end(); itc != edc; ++itc )
446  {
447 
448  TEveStraightLineSet *scposition = new TEveStraightLineSet;
449  scposition->SetDepthTest( false );
450  scposition->SetPickable( kTRUE );
451 
452  short firststrip = itc->firstStrip();
453 
454  if( &*itc == cluster )
455  {
456  scposition->SetTitle( Form( "Exact SiStripCluster from TrackingRecHit, first strip %d", firststrip ));
457  scposition->SetLineColor( kGreen );
458  }
459  else
460  {
461  scposition->SetTitle( Form( "SiStripCluster, first strip %d", firststrip ));
462  scposition->SetLineColor( kRed );
463  }
464 
465  float localTop[3] = { 0.0, 0.0, 0.0 };
466  float localBottom[3] = { 0.0, 0.0, 0.0 };
467 
468  fireworks::localSiStrip( firststrip, localTop, localBottom, pars, rawid );
469 
470  float globalTop[3];
471  float globalBottom[3];
472  geom->localToGlobal( rawid, localTop, globalTop, localBottom, globalBottom );
473 
474  scposition->AddLine( globalTop[0], globalTop[1], globalTop[2],
475  globalBottom[0], globalBottom[1], globalBottom[2] );
476 
477  setupAddElement( scposition, tList, iItem, master, false );
478  }
479  }
480  else
481  {
482  short firststrip = cluster->firstStrip();
483  TEveStraightLineSet *scposition = new TEveStraightLineSet;
484  scposition->SetDepthTest( false );
485  scposition->SetPickable( kTRUE );
486  scposition->SetTitle( Form( "SiStripCluster, first strip %d", firststrip ));
487 
488  float localTop[3] = { 0.0, 0.0, 0.0 };
489  float localBottom[3] = { 0.0, 0.0, 0.0 };
490 
491  fireworks::localSiStrip( firststrip, localTop, localBottom, pars, rawid );
492 
493  float globalTop[3];
494  float globalBottom[3];
495  geom->localToGlobal( rawid, localTop, globalTop, localBottom, globalBottom );
496 
497  scposition->AddLine( globalTop[0], globalTop[1], globalTop[2],
498  globalBottom[0], globalBottom[1], globalBottom[2] );
499 
500  setupAddElement( scposition, tList, iItem, master, true );
501  }
502  }
503  else if( !rechit->isValid() && ( rawid != 0 )) // lost hit
504  {
505  if( allClusters != 0 )
506  {
507  edmNew::DetSetVector<SiStripCluster>::const_iterator itds = allClusters->find( rawid );
508  if( itds != allClusters->end())
509  {
510  const edmNew::DetSet<SiStripCluster> & clustersOnThisDet = *itds;
511  for( edmNew::DetSet<SiStripCluster>::const_iterator itc = clustersOnThisDet.begin(), edc = clustersOnThisDet.end(); itc != edc; ++itc )
512  {
513  short firststrip = itc->firstStrip();
514 
515  TEveStraightLineSet *scposition = new TEveStraightLineSet;
516  scposition->SetDepthTest( false );
517  scposition->SetPickable( kTRUE );
518  scposition->SetTitle( Form( "Lost SiStripCluster, first strip %d", firststrip ));
519 
520  float localTop[3] = { 0.0, 0.0, 0.0 };
521  float localBottom[3] = { 0.0, 0.0, 0.0 };
522 
523  fireworks::localSiStrip( firststrip, localTop, localBottom, pars, rawid );
524 
525  float globalTop[3];
526  float globalBottom[3];
527  geom->localToGlobal( rawid, localTop, globalTop, localBottom, globalBottom );
528 
529  scposition->AddLine( globalTop[0], globalTop[1], globalTop[2],
530  globalBottom[0], globalBottom[1], globalBottom[2] );
531 
532 
533  setupAddElement( scposition, tList, iItem, master, false );
534  scposition->SetLineColor( kRed );
535  }
536  }
537  }
538  }
539  else
540  {
542  << "*ANOTHER* option possible: valid=" << rechit->isValid()
543  << ", rawid=" << rawid << std::endl;
544  }
545  }
546 }
547 
548 //______________________________________________________________________________
549 
550 void
551 pushNearbyPixelHits( std::vector<TVector3> &pixelPoints, const FWEventItem &iItem, const reco::Track &t )
552 {
553  const edmNew::DetSetVector<SiPixelCluster> * allClusters = 0;
554  for( trackingRecHit_iterator it = t.recHitsBegin(), itEnd = t.recHitsEnd(); it != itEnd; ++it)
555  {
556  if( typeid(**it) == typeid( SiPixelRecHit ))
557  {
558  const SiPixelRecHit &hit = static_cast<const SiPixelRecHit &>(**it);
559  if( hit.cluster().isNonnull() && hit.cluster().isAvailable())
560  {
562  iItem.getEvent()->get(hit.cluster().id(), allClustersHandle);
563  allClusters = allClustersHandle.product();
564  break;
565  }
566  }
567  }
568  if( allClusters == 0 ) return;
569 
570  const FWGeometry *geom = iItem.getGeom();
571 
572  for( trackingRecHit_iterator it = t.recHitsBegin(), itEnd = t.recHitsEnd(); it != itEnd; ++it )
573  {
574  const TrackingRecHit* rh = &(**it);
575 
576  DetId id = (*it)->geographicalId();
577  if( ! geom->contains( id ))
578  {
580  << "failed to get geometry of Tracker Det with raw id: "
581  << id.rawId() << std::endl;
582 
583  continue;
584  }
585 
586  // -- in which detector are we?
587  unsigned int subdet = (unsigned int)id.subdetId();
588  if(( subdet != PixelSubdetector::PixelBarrel ) && ( subdet != PixelSubdetector::PixelEndcap )) continue;
589 
590  const SiPixelCluster* hitCluster = 0;
591  if( const SiPixelRecHit* pixel = dynamic_cast<const SiPixelRecHit*>( rh ))
592  hitCluster = pixel->cluster().get();
593  edmNew::DetSetVector<SiPixelCluster>::const_iterator itds = allClusters->find(id.rawId());
594  if( itds != allClusters->end())
595  {
596  const edmNew::DetSet<SiPixelCluster> & clustersOnThisDet = *itds;
597  for( edmNew::DetSet<SiPixelCluster>::const_iterator itc = clustersOnThisDet.begin(), edc = clustersOnThisDet.end(); itc != edc; ++itc )
598  {
599  if( &*itc != hitCluster )
600  pushPixelCluster( pixelPoints, *geom, id, *itc, geom->getParameters( id ));
601  }
602  }
603  }
604 }
605 
606 //______________________________________________________________________________
607 
608 void
609 pushPixelHits( std::vector<TVector3> &pixelPoints, const FWEventItem &iItem, const reco::Track &t )
610 {
611  /*
612  * -- return for each Pixel Hit a 3D point
613  */
614  const FWGeometry *geom = iItem.getGeom();
615 
616  double dz = t.dz();
617  double vz = t.vz();
618  double etaT = t.eta();
619 
620  fwLog( fwlog::kDebug ) << "Track eta: " << etaT << ", vz: " << vz << ", dz: " << dz
621  << std::endl;
622 
623  int cnt = 0;
624  for( trackingRecHit_iterator it = t.recHitsBegin(), itEnd = t.recHitsEnd(); it != itEnd; ++it )
625  {
626  const TrackingRecHit* rh = &(**it);
627  // -- get position of center of wafer, assuming (0,0,0) is the center
628  DetId id = (*it)->geographicalId();
629  if( ! geom->contains( id ))
630  {
632  << "failed to get geometry of Tracker Det with raw id: "
633  << id.rawId() << std::endl;
634 
635  continue;
636  }
637 
638  // -- in which detector are we?
639  unsigned int subdet = (unsigned int)id.subdetId();
640 
641  if(( subdet == PixelSubdetector::PixelBarrel ) || ( subdet == PixelSubdetector::PixelEndcap ))
642  {
643  fwLog( fwlog::kDebug ) << cnt++ << " -- "
644  << subdets[subdet];
645 
646  if( const SiPixelRecHit* pixel = dynamic_cast<const SiPixelRecHit*>( rh ))
647  {
648  const SiPixelCluster& c = *( pixel->cluster());
649  pushPixelCluster( pixelPoints, *geom, id, c, geom->getParameters( id ));
650  }
651  }
652  }
653 }
654 
655 void
656 pushPixelCluster( std::vector<TVector3> &pixelPoints, const FWGeometry &geom, DetId id, const SiPixelCluster &c, const float* pars )
657 {
658  double row = c.minPixelRow();
659  double col = c.minPixelCol();
660  float lx = 0.;
661  float ly = 0.;
662 
663  int nrows = (int)pars[0];
664  int ncols = (int)pars[1];
665  lx = pixelLocalX( row, nrows );
666  ly = pixelLocalY( col, ncols );
667 
669  << ", row: " << row << ", col: " << col
670  << ", lx: " << lx << ", ly: " << ly ;
671 
672  float local[3] = { lx, ly, 0. };
673  float global[3];
674  geom.localToGlobal( id, local, global );
675  TVector3 pb( global[0], global[1], global[2] );
676  pixelPoints.push_back( pb );
677 
679  << " x: " << pb.X()
680  << ", y: " << pb.Y()
681  << " z: " << pb.Z()
682  << " eta: " << pb.Eta()
683  << ", phi: " << pb.Phi()
684  << " rho: " << pb.Pt() << std::endl;
685 }
686 
687 //______________________________________________________________________________
688 
690 info(const DetId& id) {
691  std::ostringstream oss;
692 
693  oss << "DetId: " << id.rawId() << "\n";
694 
695  switch ( id.det() ) {
696 
697  case DetId::Tracker:
698  switch ( id.subdetId() ) {
700  {
701  oss <<"TIB "<<TIBDetId(id).layer();
702  }
703  break;
705  {
706  oss <<"TOB "<<TOBDetId(id).layer();
707  }
708  break;
710  {
711  oss <<"TEC "<<TECDetId(id).wheel();
712  }
713  break;
715  {
716  oss <<"TID "<<TIDDetId(id).wheel();
717  }
718  break;
720  {
721  oss <<"PixBarrel "<< PXBDetId(id).layer();
722  }
723  break;
725  {
726  oss <<"PixEndcap "<< PXBDetId(id).layer();
727  }
728  break;
729  }
730  break;
731 
732  case DetId::Muon:
733  switch ( id.subdetId() ) {
734  case MuonSubdetId::DT:
735  {
736  DTChamberId detId(id.rawId());
737  oss << "DT chamber (wheel, station, sector): "
738  << detId.wheel() << ", "
739  << detId.station() << ", "
740  << detId.sector();
741  }
742  break;
743  case MuonSubdetId::CSC:
744  {
745  CSCDetId detId(id.rawId());
746  oss << "CSC chamber (endcap, station, ring, chamber, layer): "
747  << detId.endcap() << ", "
748  << detId.station() << ", "
749  << detId.ring() << ", "
750  << detId.chamber() << ", "
751  << detId.layer();
752  }
753  break;
754  case MuonSubdetId::RPC:
755  {
756  RPCDetId detId(id.rawId());
757  oss << "RPC chamber ";
758  switch ( detId.region() ) {
759  case 0:
760  oss << "/ barrel / (wheel, station, sector, layer, subsector, roll): "
761  << detId.ring() << ", "
762  << detId.station() << ", "
763  << detId.sector() << ", "
764  << detId.layer() << ", "
765  << detId.subsector() << ", "
766  << detId.roll();
767  break;
768  case 1:
769  oss << "/ forward endcap / (wheel, station, sector, layer, subsector, roll): "
770  << detId.ring() << ", "
771  << detId.station() << ", "
772  << detId.sector() << ", "
773  << detId.layer() << ", "
774  << detId.subsector() << ", "
775  << detId.roll();
776  break;
777  case -1:
778  oss << "/ backward endcap / (wheel, station, sector, layer, subsector, roll): "
779  << detId.ring() << ", "
780  << detId.station() << ", "
781  << detId.sector() << ", "
782  << detId.layer() << ", "
783  << detId.subsector() << ", "
784  << detId.roll();
785  break;
786  }
787  }
788  break;
789  }
790  break;
791 
792  case DetId::Calo:
793  {
794  CaloTowerDetId detId( id.rawId() );
795  oss << "CaloTower (ieta, iphi): "
796  << detId.ieta() << ", "
797  << detId.iphi();
798  }
799  break;
800 
801  case DetId::Ecal:
802  switch ( id.subdetId() ) {
803  case EcalBarrel:
804  {
805  EBDetId detId(id);
806  oss << "EcalBarrel (ieta, iphi, tower_ieta, tower_iphi): "
807  << detId.ieta() << ", "
808  << detId.iphi() << ", "
809  << detId.tower_ieta() << ", "
810  << detId.tower_iphi();
811  }
812  break;
813  case EcalEndcap:
814  {
815  EEDetId detId(id);
816  oss << "EcalEndcap (ix, iy, SuperCrystal, crystal, quadrant): "
817  << detId.ix() << ", "
818  << detId.iy() << ", "
819  << detId.isc() << ", "
820  << detId.ic() << ", "
821  << detId.iquadrant();
822  }
823  break;
824  case EcalPreshower:
825  oss << "EcalPreshower";
826  break;
827  case EcalTriggerTower:
828  oss << "EcalTriggerTower";
829  break;
830  case EcalLaserPnDiode:
831  oss << "EcalLaserPnDiode";
832  break;
833  }
834  break;
835 
836  case DetId::Hcal:
837  {
838  HcalDetId detId(id);
839  switch ( detId.subdet() ) {
840  case HcalEmpty:
841  oss << "HcalEmpty ";
842  break;
843  case HcalBarrel:
844  oss << "HcalBarrel ";
845  break;
846  case HcalEndcap:
847  oss << "HcalEndcap ";
848  break;
849  case HcalOuter:
850  oss << "HcalOuter ";
851  break;
852  case HcalForward:
853  oss << "HcalForward ";
854  break;
855  case HcalTriggerTower:
856  oss << "HcalTriggerTower ";
857  break;
858  case HcalOther:
859  oss << "HcalOther ";
860  break;
861  }
862  oss << "(ieta, iphi, depth):"
863  << detId.ieta() << ", "
864  << detId.iphi() << ", "
865  << detId.depth();
866  }
867  break;
868  default :;
869  }
870  return oss.str();
871 }
872 
874 info(const std::set<DetId>& idSet) {
876  for(std::set<DetId>::const_iterator id = idSet.begin(), idEnd = idSet.end(); id != idEnd; ++id)
877  {
878  text += info(*id);
879  text += "\n";
880  }
881  return text;
882 }
883 
885 info(const std::vector<DetId>& idSet) {
887  for(std::vector<DetId>::const_iterator id = idSet.begin(), idEnd = idSet.end(); id != idEnd; ++id)
888  {
889  text += info(*id);
890  text += "\n";
891  }
892  return text;
893 }
894 }
ClusterRef cluster() const
dbl * delta
Definition: mlp_gen.cc:36
bool isAvailable() const
Definition: Ref.h:614
int chamber() const
Definition: CSCDetId.h:81
int i
Definition: DBlmapReader.cc:9
const FWDisplayProperties & defaultDisplayProperties() const
Definition: FWEventItem.cc:451
int minPixelCol() const
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
const_iterator end(bool update=false) const
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:250
int ix() const
Definition: EEDetId.h:76
int ic() const
Definition: EEDetId.cc:324
unsigned int layer() const
layer id
Definition: TOBDetId.h:39
list parent
Definition: dbtoconf.py:74
HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:30
const TrackExtraRef & extra() const
reference to &quot;extra&quot; object
Definition: Track.h:189
CaloTopology const * topology(0)
TEveTrack * prepareTrack(const reco::Track &track, TEveTrackPropagator *propagator, const std::vector< TEveVector > &extraRefPoints=std::vector< TEveVector >())
Definition: TrackUtils.cc:73
float pixelLocalY(const double mpy, const int m_ncols)
Definition: TrackUtils.cc:201
static const int COLS_PER_ROC
Definition: TrackUtils.cc:65
int tower_ieta() const
get the HCAL/trigger ieta of this crystal
Definition: EBDetId.h:55
const float * getParameters(unsigned int id) const
Definition: FWGeometry.cc:296
bool get(ProductID const &, Handle< T > &) const
Definition: EventBase.h:104
int tower_iphi() const
get the HCAL/trigger iphi of this crystal
Definition: EBDetId.cc:114
int isc() const
Definition: EEDetId.cc:285
assert(m_qm.get())
static const int BIG_PIX_PER_ROC_Y
Definition: TrackUtils.cc:64
bool innerOk() const
return true if the innermost hit is valid
Definition: Track.h:50
int iquadrant() const
Definition: EEDetId.cc:264
uint16_t firstStrip() const
data_type const * const_iterator
Definition: DetSetNew.h:30
double px() const
x coordinate of momentum vector
Definition: TrackBase.h:675
const math::XYZPoint & outerPosition() const
position of the outermost hit
Definition: Track.h:65
int layer() const
Definition: CSCDetId.h:74
ProductID id() const
Accessor for product ID.
Definition: Ref.h:256
unsigned int layer() const
layer id
Definition: PXBDetId.h:35
Color_t color() const
static const int ROWS_PER_ROC
Definition: TrackUtils.cc:66
void localSiStrip(short strip, float *localTop, float *localBottom, const float *pars, unsigned int id)
Definition: TrackUtils.cc:295
int iphi() const
get the crystal iphi
Definition: EBDetId.h:53
const math::XYZPoint & innerPosition() const
position of the innermost hit
Definition: Track.h:55
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
int endcap() const
Definition: CSCDetId.h:106
static const double PITCHX
Definition: TrackUtils.cc:62
static const int BIG_PIX_PER_ROC_X
Definition: TrackUtils.cc:67
int depth() const
get the tower depth
Definition: HcalDetId.h:40
Char_t transparency() const
static const int CSC
Definition: MuonSubdetId.h:13
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:699
int roll() const
Definition: RPCDetId.h:120
int minPixelRow() const
int ring() const
Definition: RPCDetId.h:72
double pt() const
track transverse momentum
Definition: TrackBase.h:669
const SiStripCluster * extractClusterFromTrackingRecHit(const TrackingRecHit *rh)
Definition: TrackUtils.cc:362
static const double MICRON
Definition: TrackUtils.cc:61
int ieta() const
get the cell ieta
Definition: HcalDetId.h:36
void localToGlobal(unsigned int id, const float *local, float *global, bool translatep=true) const
Definition: FWGeometry.cc:328
int iphi() const
get the tower iphi
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
int j
Definition: DBlmapReader.cc:9
double f[11][100]
static const double PITCHY
Definition: TrackUtils.cc:63
ClusterRef cluster() const
tuple text
Definition: runonSM.py:42
int iy() const
Definition: EEDetId.h:82
bool contains(unsigned int id) const
Definition: FWGeometry.h:109
int ieta() const
get the crystal ieta
Definition: EBDetId.h:51
math::XYZPoint Point
point in the space
Definition: TrackBase.h:82
void pushPixelHits(std::vector< TVector3 > &pixelPoints, const FWEventItem &iItem, const reco::Track &t)
Definition: TrackUtils.cc:609
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:104
void addSiStripClusters(const FWEventItem *iItem, const reco::Track &t, class TEveElement *tList, bool addNearbyClusters, bool master)
Definition: TrackUtils.cc:385
double pz() const
z coordinate of momentum vector
Definition: TrackBase.h:687
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:657
int iphi() const
get the cell iphi
Definition: HcalDetId.h:38
double vz() const
z coordinate of the reference point on track
Definition: TrackBase.h:717
int ring() const
Definition: CSCDetId.h:88
const edm::EventBase * getEvent() const
Definition: FWEventItem.h:146
int layer() const
Definition: RPCDetId.h:108
Definition: DetId.h:18
const math::XYZVector & outerMomentum() const
momentum vector at the outermost hit position
Definition: Track.h:70
T const * product() const
Definition: Handle.h:81
bool outerOk() const
return true if the outermost hit is valid
Definition: Track.h:45
#define fwLog(_level_)
Definition: fwLog.h:50
ClusterRef cluster() const
Definition: SiPixelRecHit.h:49
bool isValid() const
unsigned int wheel() const
wheel id
Definition: TECDetId.h:52
int sector() const
Sector id: the group of chambers at same phi (and increasing r)
Definition: RPCDetId.h:102
unsigned int layer() const
layer id
Definition: TIBDetId.h:41
const_iterator find(id_type i, bool update=false) const
int subsector() const
SubSector id : some sectors are divided along the phi direction in subsectors (from 1 to 4 in Barrel...
Definition: RPCDetId.h:114
Pixel cluster – collection of neighboring pixels above threshold.
double vy() const
y coordinate of the reference point on track
Definition: TrackBase.h:711
iterator end()
Definition: DetSetNew.h:70
static const int RPC
Definition: MuonSubdetId.h:14
int sector() const
Definition: DTChamberId.h:61
const math::XYZVector & innerMomentum() const
momentum vector at the innermost hit position
Definition: Track.h:60
static int position[264][3]
Definition: ReadPGInfo.cc:509
std::string info(const DetId &)
Definition: TrackUtils.cc:690
int station() const
Definition: CSCDetId.h:99
int charge() const
track electric charge
Definition: TrackBase.h:615
static const int DT
Definition: MuonSubdetId.h:12
DetId geographicalId() const
int ieta() const
get the tower ieta
int station() const
Return the station number.
Definition: DTChamberId.h:51
int col
Definition: cuy.py:1008
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:45
static const std::string subdets[7]
Definition: TrackUtils.cc:70
void pushNearbyPixelHits(std::vector< TVector3 > &pixelPoints, const FWEventItem &iItem, const reco::Track &t)
Definition: TrackUtils.cc:551
math::XYZVector Vector
spatial vector
Definition: TrackBase.h:79
void pushPixelCluster(std::vector< TVector3 > &pixelPoints, const FWGeometry &geom, DetId id, const SiPixelCluster &c, const float *pars)
Definition: TrackUtils.cc:656
void setupAddElement(TEveElement *el, TEveElement *parent, const FWEventItem *item, bool master, bool color)
Definition: TrackUtils.cc:340
float pixelLocalX(const double mpx, const int m_nrows)
Definition: TrackUtils.cc:168
const FWGeometry * getGeom() const
Definition: FWEventItem.cc:681
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5
double py() const
y coordinate of momentum vector
Definition: TrackBase.h:681
TrackingRecHitCollection::base::const_iterator trackingRecHit_iterator
iterator over a vector of reference to TrackingRecHit in the same collection
double vx() const
x coordinate of the reference point on track
Definition: TrackBase.h:705
unsigned int wheel() const
wheel id
Definition: TIDDetId.h:50
Our base class.
Definition: SiPixelRecHit.h:23
trackingRecHit_iterator recHitsEnd() const
Iterator to last hit on the track.
Definition: Track.h:109
int region() const
Region id: 0 for Barrel, +/-1 For +/- Endcap.
Definition: RPCDetId.h:63
int station() const
Definition: RPCDetId.h:96
iterator begin()
Definition: DetSetNew.h:67