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 // Transform measurement to local coordinates in X dimension
260 float
261 phase2PixelLocalX( const double mpx, const float* par, const float* shape )
262 {
263  // The final position in local coordinates
264  return (-shape[1] + mpx * par[0]);
265 }
266 
267 // Transform measurement to local coordinates in Y dimension
268 float
269 phase2PixelLocalY( const double mpy, const float* par, const float* shape )
270 {
271  // The final position in local coordinates
272  return (-shape[2] + mpy * par[1]);
273 }
274 
275 //
276 // Returns strip geometry in local coordinates of a detunit.
277 // The strip is a line from a localTop to a localBottom point.
278 void localSiStrip( short strip, float* localTop, float* localBottom, const float* pars, unsigned int id )
279 {
280  Float_t topology = pars[0];
281  Float_t halfStripLength = pars[2] * 0.5;
282 
283  Double_t localCenter[3] = { 0.0, 0.0, 0.0 };
284  localTop[1] = halfStripLength;
285  localBottom[1] = -halfStripLength;
286 
287  if( topology == 1 ) // RadialStripTopology
288  {
289  // stripAngle = phiOfOneEdge + strip * angularWidth
290  // localY = originToIntersection * tan( stripAngle )
291  Float_t stripAngle = tan( pars[5] + strip * pars[6] );
292  Float_t delta = halfStripLength * stripAngle;
293  localCenter[0] = pars[4] * stripAngle;
294  localTop[0] = localCenter[0] + delta;
295  localBottom[0] = localCenter[0] - delta;
296  }
297  else if( topology == 2 ) // RectangularStripTopology
298  {
299  // offset = -numberOfStrips/2. * pitch
300  // localY = strip * pitch + offset
301  Float_t offset = -pars[1] * 0.5 * pars[3];
302  localCenter[0] = strip * pars[3] + offset;
303  localTop[0] = localCenter[0];
304  localBottom[0] = localCenter[0];
305  }
306  else if( topology == 3 ) // TrapezoidalStripTopology
307  {
309  << "did not expect TrapezoidalStripTopology of "
310  << id << std::endl;
311  }
312  else if( pars[0] == 0 ) // StripTopology
313  {
315  << "did not find StripTopology of "
316  << id << std::endl;
317  }
318 }
319 
320 //______________________________________________________________________________
321 
322 void
323 setupAddElement(TEveElement* el, TEveElement* parent, const FWEventItem* item, bool master, bool color)
324 {
325  if (master)
326  {
327  el->CSCTakeAnyParentAsMaster();
328  el->SetPickable(true);
329  }
330 
331  if (color)
332  {
333  el->CSCApplyMainColorToMatchingChildren();
334  el->CSCApplyMainTransparencyToMatchingChildren();
335  el->SetMainColor(item->defaultDisplayProperties().color());
336  assert((item->defaultDisplayProperties().transparency() >= 0)
337  && (item->defaultDisplayProperties().transparency() <= 100));
338  el->SetMainTransparency(item->defaultDisplayProperties().transparency());
339  }
340  parent->AddElement(el);
341 }
342 
343 //______________________________________________________________________________
344 
346 {
347  const SiStripCluster* cluster = nullptr;
348 
349  if( const SiStripRecHit2D* hit2D = dynamic_cast<const SiStripRecHit2D*>( rechit ))
350  {
351  fwLog( fwlog::kDebug ) << "hit 2D ";
352 
353  cluster = hit2D->cluster().get();
354  }
355  if( cluster == nullptr )
356  {
357  if( const SiStripRecHit1D* hit1D = dynamic_cast<const SiStripRecHit1D*>( rechit ))
358  {
359  fwLog( fwlog::kDebug ) << "hit 1D ";
360 
361  cluster = hit1D->cluster().get();
362  }
363  }
364  return cluster;
365 }
366 
367 void
368 addSiStripClusters( const FWEventItem* iItem, const reco::Track &t, class TEveElement *tList, bool addNearbyClusters, bool master )
369 {
370  // master is true if the product is for proxy builder
371  const FWGeometry *geom = iItem->getGeom();
372 
373  const edmNew::DetSetVector<SiStripCluster> * allClusters = nullptr;
374  if( addNearbyClusters )
375  {
376  for( trackingRecHit_iterator it = t.recHitsBegin(), itEnd = t.recHitsEnd(); it != itEnd; ++it )
377  {
378  const auto & rhs = *(*(it));
379  if( typeid( rhs ) == typeid( SiStripRecHit2D ))
380  {
381  const SiStripRecHit2D &hit = static_cast<const SiStripRecHit2D &>( **it );
382  if( hit.cluster().isNonnull() && hit.cluster().isAvailable())
383  {
385  iItem->getEvent()->get(hit.cluster().id(), allClustersHandle);
386  allClusters = allClustersHandle.product();
387  break;
388  }
389  }
390  else if( typeid( rhs ) == typeid( SiStripRecHit1D ))
391  {
392  const SiStripRecHit1D &hit = static_cast<const SiStripRecHit1D &>( **it );
393  if( hit.cluster().isNonnull() && hit.cluster().isAvailable())
394  {
396  iItem->getEvent()->get(hit.cluster().id(), allClustersHandle);
397  allClusters = allClustersHandle.product();
398  break;
399  }
400  }
401  }
402  }
403 
404  for( trackingRecHit_iterator it = t.recHitsBegin(), itEnd = t.recHitsEnd(); it != itEnd; ++it )
405  {
406  unsigned int rawid = (*it)->geographicalId();
407  if( ! geom->contains( rawid ))
408  {
410  << "failed to get geometry of SiStripCluster with detid: "
411  << rawid << std::endl;
412 
413  continue;
414  }
415 
416  const float* pars = geom->getParameters( rawid );
417 
418  // -- get phi from SiStripHit
419  auto rechitRef = *it;
420  const TrackingRecHit *rechit = &( *rechitRef );
421  const SiStripCluster *cluster = extractClusterFromTrackingRecHit( rechit );
422 
423  if( cluster )
424  {
425  if( allClusters != nullptr )
426  {
427  const edmNew::DetSet<SiStripCluster> & clustersOnThisDet = (*allClusters)[rechit->geographicalId().rawId()];
428 
429  for( edmNew::DetSet<SiStripCluster>::const_iterator itc = clustersOnThisDet.begin(), edc = clustersOnThisDet.end(); itc != edc; ++itc )
430  {
431 
432  TEveStraightLineSet *scposition = new TEveStraightLineSet;
433  scposition->SetDepthTest( false );
434  scposition->SetPickable( kTRUE );
435 
436  short firststrip = itc->firstStrip();
437 
438  if( &*itc == cluster )
439  {
440  scposition->SetTitle( Form( "Exact SiStripCluster from TrackingRecHit, first strip %d", firststrip ));
441  scposition->SetLineColor( kGreen );
442  }
443  else
444  {
445  scposition->SetTitle( Form( "SiStripCluster, first strip %d", firststrip ));
446  scposition->SetLineColor( kRed );
447  }
448 
449  float localTop[3] = { 0.0, 0.0, 0.0 };
450  float localBottom[3] = { 0.0, 0.0, 0.0 };
451 
452  fireworks::localSiStrip( firststrip, localTop, localBottom, pars, rawid );
453 
454  float globalTop[3];
455  float globalBottom[3];
456  geom->localToGlobal( rawid, localTop, globalTop, localBottom, globalBottom );
457 
458  scposition->AddLine( globalTop[0], globalTop[1], globalTop[2],
459  globalBottom[0], globalBottom[1], globalBottom[2] );
460 
461  setupAddElement( scposition, tList, iItem, master, false );
462  }
463  }
464  else
465  {
466  short firststrip = cluster->firstStrip();
467  TEveStraightLineSet *scposition = new TEveStraightLineSet;
468  scposition->SetDepthTest( false );
469  scposition->SetPickable( kTRUE );
470  scposition->SetTitle( Form( "SiStripCluster, first strip %d", firststrip ));
471 
472  float localTop[3] = { 0.0, 0.0, 0.0 };
473  float localBottom[3] = { 0.0, 0.0, 0.0 };
474 
475  fireworks::localSiStrip( firststrip, localTop, localBottom, pars, rawid );
476 
477  float globalTop[3];
478  float globalBottom[3];
479  geom->localToGlobal( rawid, localTop, globalTop, localBottom, globalBottom );
480 
481  scposition->AddLine( globalTop[0], globalTop[1], globalTop[2],
482  globalBottom[0], globalBottom[1], globalBottom[2] );
483 
484  setupAddElement( scposition, tList, iItem, master, true );
485  }
486  }
487  else if( !rechit->isValid() && ( rawid != 0 )) // lost hit
488  {
489  if( allClusters != nullptr )
490  {
491  edmNew::DetSetVector<SiStripCluster>::const_iterator itds = allClusters->find( rawid );
492  if( itds != allClusters->end())
493  {
494  const edmNew::DetSet<SiStripCluster> & clustersOnThisDet = *itds;
495  for( edmNew::DetSet<SiStripCluster>::const_iterator itc = clustersOnThisDet.begin(), edc = clustersOnThisDet.end(); itc != edc; ++itc )
496  {
497  short firststrip = itc->firstStrip();
498 
499  TEveStraightLineSet *scposition = new TEveStraightLineSet;
500  scposition->SetDepthTest( false );
501  scposition->SetPickable( kTRUE );
502  scposition->SetTitle( Form( "Lost SiStripCluster, first strip %d", firststrip ));
503 
504  float localTop[3] = { 0.0, 0.0, 0.0 };
505  float localBottom[3] = { 0.0, 0.0, 0.0 };
506 
507  fireworks::localSiStrip( firststrip, localTop, localBottom, pars, rawid );
508 
509  float globalTop[3];
510  float globalBottom[3];
511  geom->localToGlobal( rawid, localTop, globalTop, localBottom, globalBottom );
512 
513  scposition->AddLine( globalTop[0], globalTop[1], globalTop[2],
514  globalBottom[0], globalBottom[1], globalBottom[2] );
515 
516 
517  setupAddElement( scposition, tList, iItem, master, false );
518  scposition->SetLineColor( kRed );
519  }
520  }
521  }
522  }
523  else
524  {
526  << "*ANOTHER* option possible: valid=" << rechit->isValid()
527  << ", rawid=" << rawid << std::endl;
528  }
529  }
530 }
531 
532 //______________________________________________________________________________
533 
534 void
535 pushNearbyPixelHits( std::vector<TVector3> &pixelPoints, const FWEventItem &iItem, const reco::Track &t )
536 {
537  const edmNew::DetSetVector<SiPixelCluster> * allClusters = nullptr;
538  for( trackingRecHit_iterator it = t.recHitsBegin(), itEnd = t.recHitsEnd(); it != itEnd; ++it)
539  {
540  const auto & rhs = *(*(it));
541  if( typeid( rhs ) == typeid( SiPixelRecHit ))
542  {
543  const SiPixelRecHit &hit = static_cast<const SiPixelRecHit &>(**it);
544  if( hit.cluster().isNonnull() && hit.cluster().isAvailable())
545  {
547  iItem.getEvent()->get(hit.cluster().id(), allClustersHandle);
548  allClusters = allClustersHandle.product();
549  break;
550  }
551  }
552  }
553  if( allClusters == nullptr ) return;
554 
555  const FWGeometry *geom = iItem.getGeom();
556 
557  for( trackingRecHit_iterator it = t.recHitsBegin(), itEnd = t.recHitsEnd(); it != itEnd; ++it )
558  {
559  const TrackingRecHit* rh = &(**it);
560 
561  DetId id = (*it)->geographicalId();
562  if( ! geom->contains( id ))
563  {
565  << "failed to get geometry of Tracker Det with raw id: "
566  << id.rawId() << std::endl;
567 
568  continue;
569  }
570 
571  // -- in which detector are we?
572  unsigned int subdet = (unsigned int)id.subdetId();
573  if(( subdet != PixelSubdetector::PixelBarrel ) && ( subdet != PixelSubdetector::PixelEndcap )) continue;
574 
575  const SiPixelCluster* hitCluster = nullptr;
576  if( const SiPixelRecHit* pixel = dynamic_cast<const SiPixelRecHit*>( rh ))
577  hitCluster = pixel->cluster().get();
578  edmNew::DetSetVector<SiPixelCluster>::const_iterator itds = allClusters->find(id.rawId());
579  if( itds != allClusters->end())
580  {
581  const edmNew::DetSet<SiPixelCluster> & clustersOnThisDet = *itds;
582  for( edmNew::DetSet<SiPixelCluster>::const_iterator itc = clustersOnThisDet.begin(), edc = clustersOnThisDet.end(); itc != edc; ++itc )
583  {
584  if( &*itc != hitCluster )
585  pushPixelCluster( pixelPoints, *geom, id, *itc, geom->getParameters( id ));
586  }
587  }
588  }
589 }
590 
591 //______________________________________________________________________________
592 
593 void
594 pushPixelHits( std::vector<TVector3> &pixelPoints, const FWEventItem &iItem, const reco::Track &t )
595 {
596  /*
597  * -- return for each Pixel Hit a 3D point
598  */
599  const FWGeometry *geom = iItem.getGeom();
600 
601  double dz = t.dz();
602  double vz = t.vz();
603  double etaT = t.eta();
604 
605  fwLog( fwlog::kDebug ) << "Track eta: " << etaT << ", vz: " << vz << ", dz: " << dz
606  << std::endl;
607 
608  int cnt = 0;
609  for( trackingRecHit_iterator it = t.recHitsBegin(), itEnd = t.recHitsEnd(); it != itEnd; ++it )
610  {
611  const TrackingRecHit* rh = &(**it);
612  // -- get position of center of wafer, assuming (0,0,0) is the center
613  DetId id = (*it)->geographicalId();
614  if( ! geom->contains( id ))
615  {
617  << "failed to get geometry of Tracker Det with raw id: "
618  << id.rawId() << std::endl;
619 
620  continue;
621  }
622 
623  // -- in which detector are we?
624  unsigned int subdet = (unsigned int)id.subdetId();
625 
626  if(( subdet == PixelSubdetector::PixelBarrel ) || ( subdet == PixelSubdetector::PixelEndcap ))
627  {
628  fwLog( fwlog::kDebug ) << cnt++ << " -- "
629  << subdets[subdet];
630 
631  if( const SiPixelRecHit* pixel = dynamic_cast<const SiPixelRecHit*>( rh ))
632  {
633  const SiPixelCluster& c = *( pixel->cluster());
634  pushPixelCluster( pixelPoints, *geom, id, c, geom->getParameters( id ));
635  }
636  }
637  }
638 }
639 
640 void
641 pushPixelCluster( std::vector<TVector3> &pixelPoints, const FWGeometry &geom, DetId id, const SiPixelCluster &c, const float* pars )
642 {
643  double row = c.minPixelRow();
644  double col = c.minPixelCol();
645  float lx = 0.;
646  float ly = 0.;
647 
648  // int nrows = (int)pars[0];
649  // int ncols = (int)pars[1];
650  lx = pixelLocalX( row, pars );
651  ly = pixelLocalY( col, pars );
652 
654  << ", row: " << row << ", col: " << col
655  << ", lx: " << lx << ", ly: " << ly ;
656 
657  float local[3] = { lx, ly, 0. };
658  float global[3];
659  geom.localToGlobal( id, local, global );
660  TVector3 pb( global[0], global[1], global[2] );
661  pixelPoints.push_back( pb );
662 
664  << " x: " << pb.X()
665  << ", y: " << pb.Y()
666  << " z: " << pb.Z()
667  << " eta: " << pb.Eta()
668  << ", phi: " << pb.Phi()
669  << " rho: " << pb.Pt() << std::endl;
670 }
671 
672 //______________________________________________________________________________
673 
675 info(const DetId& id) {
676  std::ostringstream oss;
677 
678  oss << "DetId: " << id.rawId() << "\n";
679 
680  switch ( id.det() ) {
681 
682  case DetId::Tracker:
684  break;
685 
686  case DetId::Muon:
687  switch ( id.subdetId() ) {
688  case MuonSubdetId::DT:
689  {
690  DTChamberId detId(id.rawId());
691  oss << "DT chamber (wheel, station, sector): "
692  << detId.wheel() << ", "
693  << detId.station() << ", "
694  << detId.sector();
695  }
696  break;
697  case MuonSubdetId::CSC:
698  {
699  CSCDetId detId(id.rawId());
700  oss << "CSC chamber (endcap, station, ring, chamber, layer): "
701  << detId.endcap() << ", "
702  << detId.station() << ", "
703  << detId.ring() << ", "
704  << detId.chamber() << ", "
705  << detId.layer();
706  }
707  break;
708  case MuonSubdetId::RPC:
709  {
710  RPCDetId detId(id.rawId());
711  oss << "RPC chamber ";
712  switch ( detId.region() ) {
713  case 0:
714  oss << "/ barrel / (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  case 1:
723  oss << "/ forward endcap / (wheel, station, sector, layer, subsector, roll): "
724  << detId.ring() << ", "
725  << detId.station() << ", "
726  << detId.sector() << ", "
727  << detId.layer() << ", "
728  << detId.subsector() << ", "
729  << detId.roll();
730  break;
731  case -1:
732  oss << "/ backward endcap / (wheel, station, sector, layer, subsector, roll): "
733  << detId.ring() << ", "
734  << detId.station() << ", "
735  << detId.sector() << ", "
736  << detId.layer() << ", "
737  << detId.subsector() << ", "
738  << detId.roll();
739  break;
740  }
741  }
742  break;
743  case MuonSubdetId::GEM:
744  {
745  GEMDetId detId(id.rawId());
746  oss << "GEM chamber (region, station, ring, chamber, layer): "
747  << detId.region() << ", "
748  << detId.station() << ", "
749  << detId.ring() << ", "
750  << detId.chamber() << ", "
751  << detId.layer();
752  }
753  break;
754  case MuonSubdetId::ME0:
755  {
756  ME0DetId detId(id.rawId());
757  oss << "ME0 chamber (region, chamber, layer): "
758  << detId.region() << ", "
759  << detId.chamber() << ", "
760  << detId.layer();
761  }
762  break;
763  }
764  break;
765 
766  case DetId::Calo:
767  {
768  CaloTowerDetId detId( id.rawId() );
769  oss << "CaloTower (ieta, iphi): "
770  << detId.ieta() << ", "
771  << detId.iphi();
772  }
773  break;
774 
775  case DetId::Ecal:
776  switch ( id.subdetId() ) {
777  case EcalBarrel:
778  {
779  EBDetId detId(id);
780  oss << "EcalBarrel (ieta, iphi, tower_ieta, tower_iphi): "
781  << detId.ieta() << ", "
782  << detId.iphi() << ", "
783  << detId.tower_ieta() << ", "
784  << detId.tower_iphi();
785  }
786  break;
787  case EcalEndcap:
788  {
789  EEDetId detId(id);
790  oss << "EcalEndcap (ix, iy, SuperCrystal, crystal, quadrant): "
791  << detId.ix() << ", "
792  << detId.iy() << ", "
793  << detId.isc() << ", "
794  << detId.ic() << ", "
795  << detId.iquadrant();
796  }
797  break;
798  case EcalPreshower:
799  oss << "EcalPreshower";
800  break;
801  case EcalTriggerTower:
802  oss << "EcalTriggerTower";
803  break;
804  case EcalLaserPnDiode:
805  oss << "EcalLaserPnDiode";
806  break;
807  }
808  break;
809 
810  case DetId::Hcal:
811  {
812  HcalDetId detId(id);
813  switch ( detId.subdet() ) {
814  case HcalEmpty:
815  oss << "HcalEmpty ";
816  break;
817  case HcalBarrel:
818  oss << "HcalBarrel ";
819  break;
820  case HcalEndcap:
821  oss << "HcalEndcap ";
822  break;
823  case HcalOuter:
824  oss << "HcalOuter ";
825  break;
826  case HcalForward:
827  oss << "HcalForward ";
828  break;
829  case HcalTriggerTower:
830  oss << "HcalTriggerTower ";
831  break;
832  case HcalOther:
833  oss << "HcalOther ";
834  break;
835  }
836  oss << "(ieta, iphi, depth):"
837  << detId.ieta() << ", "
838  << detId.iphi() << ", "
839  << detId.depth();
840  }
841  break;
842  default :;
843  }
844  return oss.str();
845 }
846 
848 info(const std::set<DetId>& idSet) {
850  for(std::set<DetId>::const_iterator id = idSet.begin(), idEnd = idSet.end(); id != idEnd; ++id)
851  {
852  text += info(*id);
853  text += "\n";
854  }
855  return text;
856 }
857 
859 info(const std::vector<DetId>& idSet) {
861  for(std::vector<DetId>::const_iterator id = idSet.begin(), idEnd = idSet.end(); id != idEnd; ++id)
862  {
863  text += info(*id);
864  text += "\n";
865  }
866  return text;
867 }
868 }
ClusterRef cluster() const
dbl * delta
Definition: mlp_gen.cc:36
bool isAvailable() const
Definition: Ref.h:575
int chamber() const
Definition: CSCDetId.h:68
Master< F > master(const F &f)
Definition: FunctClone.h:68
const FWDisplayProperties & defaultDisplayProperties() const
Definition: FWEventItem.cc:456
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:251
int ix() const
Definition: EEDetId.h:77
int ic() const
Definition: EEDetId.cc:245
HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:146
const TrackExtraRef & extra() const
reference to "extra" object
Definition: Track.h:194
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:53
const float * getParameters(unsigned int id) const
Definition: FWGeometry.cc:446
bool get(ProductID const &, Handle< T > &) const
Definition: EventBase.h:103
const FWGeometry * getGeom() const
Definition: Context.h:83
int tower_iphi() const
get the HCAL/trigger iphi of this crystal
Definition: EBDetId.cc:100
int isc() const
Definition: EEDetId.cc:222
static const int GEM
Definition: MuonSubdetId.h:15
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:50
int ring() const
Definition: GEMDetId.h:59
float phase2PixelLocalY(const double mpy, const float *, const float *)
Definition: TrackUtils.cc:269
bool innerOk() const
return true if the innermost hit is valid
Definition: Track.h:52
int iquadrant() const
Definition: EEDetId.cc:206
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:660
const math::XYZPoint & outerPosition() const
position of the outermost hit
Definition: Track.h:67
std::string print(DetId detid) const
int layer() const
Definition: CSCDetId.h:61
ProductID id() const
Accessor for product ID.
Definition: Ref.h:257
Color_t color() const
void localSiStrip(short strip, float *localTop, float *localBottom, const float *pars, unsigned int id)
Definition: TrackUtils.cc:278
int iphi() const
get the crystal iphi
Definition: EBDetId.h:51
const math::XYZPoint & innerPosition() const
position of the innermost hit
Definition: Track.h:57
int endcap() const
Definition: CSCDetId.h:93
float phase2PixelLocalX(const double mpx, const float *, const float *)
Definition: TrackUtils.cc:261
static const int ME0
Definition: MuonSubdetId.h:16
int depth() const
get the tower depth
Definition: HcalDetId.h:166
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:684
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:654
const SiStripCluster * extractClusterFromTrackingRecHit(const TrackingRecHit *rh)
Definition: TrackUtils.cc:345
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:159
void localToGlobal(unsigned int id, const float *local, float *global, bool translatep=true) const
Definition: FWGeometry.cc:478
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:83
bool contains(unsigned int id) const
Definition: FWGeometry.h:117
#define end
Definition: vmac.h:39
int ieta() const
get the crystal ieta
Definition: EBDetId.h:49
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:594
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:109
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:368
double pz() const
z coordinate of momentum vector
Definition: TrackBase.h:672
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:642
int iphi() const
get the cell iphi
Definition: HcalDetId.h:161
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:702
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:72
T const * product() const
Definition: Handle.h:74
bool outerOk() const
return true if the outermost hit is valid
Definition: Track.h:47
#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:696
#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:62
static int position[264][3]
Definition: ReadPGInfo.cc:509
std::string info(const DetId &)
Definition: TrackUtils.cc:675
col
Definition: cuy.py:1010
int station() const
Definition: CSCDetId.h:86
int charge() const
track electric charge
Definition: TrackBase.h:600
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:133
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:535
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:641
void setupAddElement(TEveElement *el, TEveElement *parent, const FWEventItem *item, bool master, bool color)
Definition: TrackUtils.cc:323
const FWGeometry * getGeom() const
Definition: FWEventItem.cc:686
#define constexpr
*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:666
double vx() const
x coordinate of the reference point on track
Definition: TrackBase.h:690
Our base class.
Definition: SiPixelRecHit.h:23
trackingRecHit_iterator recHitsEnd() const
Iterator to last hit on the track.
Definition: Track.h:114
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