CMS 3D CMS Logo

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