CMS 3D CMS Logo

TrackerTopology.cc
Go to the documentation of this file.
6 #include <sstream>
7 
9  const TECValues& tecv, const TIBValues& tibv,
10  const TIDValues& tidv, const TOBValues& tobv)
11  : pbVals_(pxb),
12  pfVals_(pxf),
13  tobVals_(tobv),
14  tibVals_(tibv),
15  tidVals_(tidv),
16  tecVals_(tecv),
17  bits_per_field{
26  }
27 {}
28 
29 
30 
31 unsigned int TrackerTopology::side(const DetId &id) const {
32  uint32_t subdet=id.subdetId();
33  if ( subdet == PixelSubdetector::PixelBarrel )
34  return 0;
35  if ( subdet == PixelSubdetector::PixelEndcap )
36  return pxfSide(id);
37  if ( subdet == StripSubdetector::TIB )
38  return 0;
39  if ( subdet == StripSubdetector::TID )
40  return tidSide(id);
41  if ( subdet == StripSubdetector::TOB )
42  return 0;
43  if ( subdet == StripSubdetector::TEC )
44  return tecSide(id);
45 
46  throw cms::Exception("Invalid DetId") << "Unsupported DetId in TrackerTopology::side";
47  return 0;
48 }
49 
50 unsigned int TrackerTopology::layer(const DetId &id) const {
51  uint32_t subdet=id.subdetId();
52  if ( subdet == PixelSubdetector::PixelBarrel )
53  return pxbLayer(id);
54  if ( subdet == PixelSubdetector::PixelEndcap )
55  return pxfDisk(id);
56  if ( subdet == StripSubdetector::TIB )
57  return tibLayer(id);
58  if ( subdet == StripSubdetector::TID )
59  return tidWheel(id);
60  if ( subdet == StripSubdetector::TOB )
61  return tobLayer(id);
62  if ( subdet == StripSubdetector::TEC )
63  return tecWheel(id);
64 
65  throw cms::Exception("Invalid DetId") << "Unsupported DetId in TrackerTopology::layer";
66  return 0;
67 }
68 
69 unsigned int TrackerTopology::module(const DetId &id) const {
70  uint32_t subdet=id.subdetId();
71  if ( subdet == PixelSubdetector::PixelBarrel )
72  return pxbModule(id);
73  if ( subdet == PixelSubdetector::PixelEndcap )
74  return pxfModule(id);
75  if ( subdet == StripSubdetector::TIB )
76  return tibModule(id);
77  if ( subdet == StripSubdetector::TID )
78  return tidModule(id);
79  if ( subdet == StripSubdetector::TOB )
80  return tobModule(id);
81  if ( subdet == StripSubdetector::TEC )
82  return tecModule(id);
83 
84  throw cms::Exception("Invalid DetId") << "Unsupported DetId in TrackerTopology::module";
85  return 0;
86 }
87 
88 uint32_t TrackerTopology::glued(const DetId &id) const {
89 
90  uint32_t subdet=id.subdetId();
91  if ( subdet == PixelSubdetector::PixelBarrel )
92  return 0;
93  if ( subdet == PixelSubdetector::PixelEndcap )
94  return 0;
95  if ( subdet == StripSubdetector::TIB )
96  return tibGlued(id);
97  if ( subdet == StripSubdetector::TID )
98  return tidGlued(id);
99  if ( subdet == StripSubdetector::TOB )
100  return tobGlued(id);
101  if ( subdet == StripSubdetector::TEC )
102  return tecGlued(id);
103 
104  throw cms::Exception("Invalid DetId") << "Unsupported DetId in TrackerTopology::glued";
105  return 0;
106 }
107 
108 uint32_t TrackerTopology::stack(const DetId &id) const {
109 
110  uint32_t subdet=id.subdetId();
111  if ( subdet == PixelSubdetector::PixelBarrel )
112  return 0;
113  if ( subdet == PixelSubdetector::PixelEndcap )
114  return 0;
115  if ( subdet == StripSubdetector::TIB )
116  return tibStack(id);
117  if ( subdet == StripSubdetector::TID )
118  return tidStack(id);
119  if ( subdet == StripSubdetector::TOB )
120  return tobStack(id);
121  if ( subdet == StripSubdetector::TEC )
122  return tecStack(id);
123 
124  throw cms::Exception("Invalid DetId") << "Unsupported DetId in TrackerTopology::stack";
125 }
126 
127 uint32_t TrackerTopology::lower(const DetId &id) const {
128 
129  uint32_t subdet=id.subdetId();
130  if ( subdet == PixelSubdetector::PixelBarrel )
131  return 0;
132  if ( subdet == PixelSubdetector::PixelEndcap )
133  return 0;
134  if ( subdet == StripSubdetector::TIB )
135  return tibLower(id);
136  if ( subdet == StripSubdetector::TID )
137  return tidLower(id);
138  if ( subdet == StripSubdetector::TOB )
139  return tobLower(id);
140  if ( subdet == StripSubdetector::TEC )
141  return tecLower(id);
142 
143  throw cms::Exception("Invalid DetId") << "Unsupported DetId in TrackerTopology::lower";
144 }
145 
146 uint32_t TrackerTopology::upper(const DetId &id) const {
147 
148  uint32_t subdet=id.subdetId();
149  if ( subdet == PixelSubdetector::PixelBarrel )
150  return 0;
151  if ( subdet == PixelSubdetector::PixelEndcap )
152  return 0;
153  if ( subdet == StripSubdetector::TIB )
154  return tibUpper(id);
155  if ( subdet == StripSubdetector::TID )
156  return tidUpper(id);
157  if ( subdet == StripSubdetector::TOB )
158  return tobUpper(id);
159  if ( subdet == StripSubdetector::TEC )
160  return tecUpper(id);
161 
162  throw cms::Exception("Invalid DetId") << "Unsupported DetId in TrackerTopology::upper";
163 }
164 
165 
166 bool TrackerTopology::isStereo(const DetId &id) const {
167 
168  uint32_t subdet=id.subdetId();
169  if ( subdet == PixelSubdetector::PixelBarrel )
170  return false;
171  if ( subdet == PixelSubdetector::PixelEndcap )
172  return false;
173  if ( subdet == StripSubdetector::TIB )
174  return tibStereo(id)!=0;
175  if ( subdet == StripSubdetector::TID )
176  return tidStereo(id)!=0;
177  if ( subdet == StripSubdetector::TOB )
178  return tobStereo(id)!=0;
179  if ( subdet == StripSubdetector::TEC )
180  return tecStereo(id)!=0;
181 
182  throw cms::Exception("Invalid DetId") << "Unsupported DetId in TrackerTopology::isStereo";
183  return false;
184 }
185 
186 bool TrackerTopology::isRPhi(const DetId &id) const {
187 
188  uint32_t subdet=id.subdetId();
189  if ( subdet == PixelSubdetector::PixelBarrel )
190  return false;
191  if ( subdet == PixelSubdetector::PixelEndcap )
192  return false;
193  if ( subdet == StripSubdetector::TIB )
194  return tibRPhi(id)!=0;
195  if ( subdet == StripSubdetector::TID )
196  return tidRPhi(id)!=0;
197  if ( subdet == StripSubdetector::TOB )
198  return tobRPhi(id)!=0;
199  if ( subdet == StripSubdetector::TEC )
200  return tecRPhi(id)!=0;
201 
202  throw cms::Exception("Invalid DetId") << "Unsupported DetId in TrackerTopology::isRPhi";
203  return false;
204 }
205 bool TrackerTopology::isLower(const DetId &id) const {
206 
207  uint32_t subdet=id.subdetId();
208  if ( subdet == PixelSubdetector::PixelBarrel )
209  return false;
210  if ( subdet == PixelSubdetector::PixelEndcap )
211  return false;
212  if ( subdet == StripSubdetector::TIB )
213  return tibLower(id)!=0;
214  if ( subdet == StripSubdetector::TID )
215  return tidLower(id)!=0;
216  if ( subdet == StripSubdetector::TOB )
217  return tobLower(id)!=0;
218  if ( subdet == StripSubdetector::TEC )
219  return tecLower(id)!=0;
220 
221  throw cms::Exception("Invalid DetId") << "Unsupported DetId in TrackerTopology::isLower";
222  return false;
223 
224 }
225 
226 bool TrackerTopology::isUpper(const DetId &id) const {
227 
228  uint32_t subdet=id.subdetId();
229  if ( subdet == PixelSubdetector::PixelBarrel )
230  return false;
231  if ( subdet == PixelSubdetector::PixelEndcap )
232  return false;
233  if ( subdet == StripSubdetector::TIB )
234  return tibUpper(id)!=0;
235  if ( subdet == StripSubdetector::TID )
236  return tidUpper(id)!=0;
237  if ( subdet == StripSubdetector::TOB )
238  return tobUpper(id)!=0;
239  if ( subdet == StripSubdetector::TEC )
240  return tecUpper(id)!=0;
241 
242  throw cms::Exception("Invalid DetId") << "Unsupported DetId in TrackerTopology::isUpper";
243  return false;
244 }
245 
247 
248  uint32_t subdet=id.subdetId();
249  if ( subdet == PixelSubdetector::PixelBarrel )
250  return 0;
251  if ( subdet == PixelSubdetector::PixelEndcap )
252  return 0;
253  if ( subdet == StripSubdetector::TIB )
254  return tibPartnerDetId(id);
255  if ( subdet == StripSubdetector::TID )
256  return tidPartnerDetId(id);
257  if ( subdet == StripSubdetector::TOB )
258  return tobPartnerDetId(id);
259  if ( subdet == StripSubdetector::TEC )
260  return tecPartnerDetId(id);
261 
262  throw cms::Exception("Invalid DetId") << "Unsupported DetId in TrackerTopology::partnerDetId";
263  return 0;
264 }
265 
267  uint32_t subdet=id.subdetId();
268  std::stringstream strstr;
269 
270  if ( subdet == PixelSubdetector::PixelBarrel ) {
271  unsigned int theLayer = pxbLayer(id);
272  unsigned int theLadder = pxbLadder(id);
273  unsigned int theModule = pxbModule(id);
274  strstr << "PixelBarrel"
275  << " Layer " << theLayer
276  << " Ladder " << theLadder
277  << " Module " << theModule ;
278  strstr << " (" << id.rawId() << ")";
279  return strstr.str();
280  }
281 
282  if ( subdet == PixelSubdetector::PixelEndcap ) {
283  unsigned int theSide = pxfSide(id);
284  unsigned int theDisk = pxfDisk(id);
285  unsigned int theBlade = pxfBlade(id);
286  unsigned int thePanel = pxfPanel(id);
287  unsigned int theModule = pxfModule(id);
288  std::string side = (pxfSide(id) == 1 ) ? "-" : "+";
289  strstr << "PixelEndcap"
290  << " Side " << theSide << side
291  << " Disk " << theDisk
292  << " Blade " << theBlade
293  << " Panel " << thePanel
294  << " Module " << theModule ;
295  strstr << " (" << id.rawId() << ")";
296  return strstr.str();
297  }
298 
299  if ( subdet == StripSubdetector::TIB ) {
300  unsigned int theLayer = tibLayer(id);
301  std::vector<unsigned int> theString = tibStringInfo(id);
302  unsigned int theModule = tibModule(id);
305  side = (theString[0] == 1 ) ? "-" : "+";
306  part = (theString[1] == 1 ) ? "int" : "ext";
308  type = (isStereo(id)) ? "stereo" : type;
309  type = (isRPhi(id)) ? "r-phi" : type;
310  type = (isStereo(id) || isRPhi(id)) ? type+" glued": "module";
311  std::string typeUpgrade;
312  typeUpgrade = (isLower(id)) ? "lower" : typeUpgrade;
313  typeUpgrade = (isUpper(id)) ? "upper" : typeUpgrade;
314  typeUpgrade = (isUpper(id) || isLower(id)) ? typeUpgrade+" stack": "module";
315  strstr << "TIB" << side
316  << " Layer " << theLayer << " " << part
317  << " String " << theString[2];
318  strstr << " Module for phase0 " << theModule << " " << type;
319  strstr << " Module for phase2 " << theModule << " " << typeUpgrade;
320  strstr << " (" << id.rawId() << ")";
321  return strstr.str();
322  }
323 
324  if ( subdet == StripSubdetector::TID ) {
325  unsigned int theSide = tidSide(id);
326  unsigned int theWheel = tidWheel(id);
327  unsigned int theRing = tidRing(id);
328  std::vector<unsigned int> theModule = tidModuleInfo(id);
331  side = (tidSide(id) == 1 ) ? "-" : "+";
332  part = (theModule[0] == 1 ) ? "back" : "front";
334  type = (isStereo(id)) ? "stereo" : type;
335  type = (isRPhi(id)) ? "r-phi" : type;
336  type = (isStereo(id) || isRPhi(id)) ? type+" glued": "module";
337  std::string typeUpgrade;
338  typeUpgrade = (isLower(id)) ? "lower" : typeUpgrade;
339  typeUpgrade = (isUpper(id)) ? "upper" : typeUpgrade;
340  typeUpgrade = (isUpper(id) || isLower(id)) ? typeUpgrade+" stack": "module";
341  strstr << "TID"
342  << " Side " << theSide << side
343  << " Wheel " << theWheel
344  << " Ring " << theRing << " " << part;
345  strstr << " Module for phase0 " << theModule[1] << " " << type;
346  strstr << " Module for phase2 " << theModule[1] << " " << typeUpgrade;
347  strstr << " (" << id.rawId() << ")";
348  return strstr.str();
349  }
350 
351  if ( subdet == StripSubdetector::TOB ) {
352  unsigned int theLayer = tobLayer(id);
353  std::vector<unsigned int> theRod = tobRodInfo(id);
354  unsigned int theModule = tobModule(id);
357  side = (((theRod[0] == 1 ) ? "-" : ((theRod[0] == 2 ) ? "+" : (theRod[0] == 3 ) ? "0" : "")));
358 // side = (theRod[0] == 2 ) ? "+" : "";
359 // side = (theRod[0] == 3 ) ? "0" : "";
361  type = (isStereo(id)) ? "stereo" : type;
362  type = (isRPhi(id)) ? "r-phi" : type;
363  type = (isStereo(id) || isRPhi(id)) ? type+" glued": "module";
364  std::string typeUpgrade;
365  typeUpgrade = (isLower(id)) ? "lower" : typeUpgrade;
366  typeUpgrade = (isUpper(id)) ? "upper" : typeUpgrade;
367  typeUpgrade = (isUpper(id) || isLower(id)) ? typeUpgrade+" stack": "module";
368  strstr << "TOB" << side
369  << " Layer " << theLayer
370  << " Rod " << theRod[1];
371  strstr << " Module for phase0 " << theModule << " " << type;
372  strstr << " Module for phase2 " << theModule << " " << typeUpgrade;
373  strstr << " (" << id.rawId() << ")";
374  return strstr.str();
375  }
376 
377  if ( subdet == StripSubdetector::TEC ) {
378  unsigned int theSide = tecSide(id);
379  unsigned int theWheel = tecWheel(id);
380  unsigned int theModule = tecModule(id);
381  std::vector<unsigned int> thePetal = tecPetalInfo(id);
382  unsigned int theRing = tecRing(id);
384  std::string petal;
385  side = (tecSide(id) == 1 ) ? "-" : "+";
386  petal = (thePetal[0] == 1 ) ? "back" : "front";
388  type = (isStereo(id)) ? "stereo" : type;
389  type = (isRPhi(id)) ? "r-phi" : type;
390  type = (isStereo(id) || isRPhi(id)) ? type+" glued": "module";
391  std::string typeUpgrade;
392  typeUpgrade = (isLower(id)) ? "lower" : typeUpgrade;
393  typeUpgrade = (isUpper(id)) ? "upper" : typeUpgrade;
394  typeUpgrade = (isUpper(id) || isLower(id)) ? typeUpgrade+" stack": "module";
395  strstr << "TEC"
396  << " Side " << theSide << side
397  << " Wheel " << theWheel
398  << " Petal " << thePetal[1] << " " << petal
399  << " Ring " << theRing;
400  strstr << " Module for phase0 " << theModule << " " << type;
401  strstr << " Module for phase2 " << theModule << " " << typeUpgrade;
402  strstr << " (" << id.rawId() << ")";
403 
404  return strstr.str();
405  }
406 
407 
408  throw cms::Exception("Invalid DetId") << "Unsupported DetId in TrackerTopology::module";
409  return strstr.str();
410 }
411 
412 
414  switch(id.subdetId()) {
417  case StripSubdetector::TID: switch (tidRing(id)) {
418  case 1: return SiStripDetId::W1A;
419  case 2: return SiStripDetId::W2A;
420  case 3: return SiStripDetId::W3A;
421  }
422  case StripSubdetector::TEC: switch (tecRing(id)) {
423  case 1: return SiStripDetId::W1B;
424  case 2: return SiStripDetId::W2B;
425  case 3: return SiStripDetId::W3B;
426  case 4: return SiStripDetId::W4;
427  //generic function to return DetIds and boolean factors
428  case 5: return SiStripDetId::W5;
429  case 6: return SiStripDetId::W6;
430  case 7: return SiStripDetId::W7;
431  }
432  }
434 }
436  int layer = -1;
437 
438  if (id.det() == DetId::Tracker) {
439  if (id.subdetId() == StripSubdetector::TOB) {
440  layer = tobLayer(id);
441  } else if (id.subdetId() == StripSubdetector::TID) {
442  layer = 100 * tidSide(id) + tidWheel(id);
443  } else {
444  edm::LogInfo("TrackerTopology") << ">>> Invalid subdetId() " ;
445  }
446  }
447  return layer;
448 }
449 
451  int layer = -1;
452 
453  if (id.det() == DetId::Tracker) {
454  if (id.subdetId() == PixelSubdetector::PixelBarrel) {
455  layer = pxbLayer(id);
456  } else if (id.subdetId() == PixelSubdetector::PixelEndcap) {
457  layer = 100 * pxfSide(id) + pxfDisk(id);
458  } else {
459  edm::LogInfo("TrackerTopology") << ">>> Invalid subdetId() " ;
460  }
461  }
462  return layer;
463 }
464 
DetId tidPartnerDetId(const DetId &id) const
type
Definition: HCALResponse.h:21
uint32_t tibStack(const DetId &id) const
uint32_t tobGlued(const DetId &id) const
uint32_t upper(const DetId &id) const
unsigned int tibLayer(const DetId &id) const
unsigned int tidRing(const DetId &id) const
uint32_t tobLower(const DetId &id) const
const PixelEndcapValues pfVals_
std::vector< unsigned int > tidModuleInfo(const DetId &id) const
unsigned int pxfDisk(const DetId &id) const
unsigned int tecRing(const DetId &id) const
ring id
uint32_t tobStereo(const DetId &id) const
unsigned int pxbLadder(const DetId &id) const
uint32_t tidRPhi(const DetId &id) const
uint32_t tobStack(const DetId &id) const
unsigned int side(const DetId &id) const
uint32_t tecGlued(const DetId &id) const
bool isStereo(const DetId &id) const
unsigned int tidWheel(const DetId &id) const
unsigned int pxbModule(const DetId &id) const
SiStripDetId::ModuleGeometry moduleGeometry(const DetId &id) const
uint32_t tecRPhi(const DetId &id) const
std::vector< unsigned int > tibStringInfo(const DetId &id) const
std::string print(DetId detid) const
uint32_t tecPartnerDetId(const DetId &id) const
DetId partnerDetId(const DetId &id) const
DetId tibPartnerDetId(const DetId &id) const
DetId tobPartnerDetId(const DetId &id) const
const PixelBarrelValues pbVals_
bool isLower(const DetId &id) const
unsigned int module(const DetId &id) const
uint32_t tibRPhi(const DetId &id) const
std::vector< unsigned int > tecPetalInfo(const DetId &id) const
uint32_t tecUpper(const DetId &id) const
unsigned int tidSide(const DetId &id) const
uint32_t tidStereo(const DetId &id) const
std::vector< unsigned int > tobRodInfo(const DetId &id) const
unsigned int tidModule(const DetId &id) const
uint32_t tidStack(const DetId &id) const
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
int getITPixelLayerNumber(const DetId &id) const
unsigned int tibModule(const DetId &id) const
unsigned int pxfModule(const DetId &id) const
unsigned int pxbLayer(const DetId &id) const
bool isUpper(const DetId &id) const
unsigned int tecModule(const DetId &id) const
Definition: DetId.h:18
uint32_t tibGlued(const DetId &id) const
uint32_t stack(const DetId &id) const
part
Definition: HCALResponse.h:20
uint32_t tidLower(const DetId &id) const
uint32_t tecLower(const DetId &id) const
unsigned int tobModule(const DetId &id) const
unsigned int layer(const DetId &id) const
uint32_t lower(const DetId &id) const
uint32_t tecStereo(const DetId &id) const
uint32_t tibUpper(const DetId &id) const
uint32_t tecStack(const DetId &id) const
unsigned int pxfSide(const DetId &id) const
uint32_t tibLower(const DetId &id) const
uint32_t tobUpper(const DetId &id) const
bool isRPhi(const DetId &id) const
uint32_t tibStereo(const DetId &id) const
uint32_t glued(const DetId &id) const
unsigned int tecWheel(const DetId &id) const
uint32_t tidUpper(const DetId &id) const
TrackerTopology(const PixelBarrelValues &pxb, const PixelEndcapValues &pxf, const TECValues &tecv, const TIBValues &tibv, const TIDValues &tidv, const TOBValues &tobv)
unsigned int pxfPanel(const DetId &id) const
unsigned int pxfBlade(const DetId &id) const
int getOTLayerNumber(const DetId &id) const
unsigned int tobLayer(const DetId &id) const
unsigned int tecSide(const DetId &id) const
uint32_t tobRPhi(const DetId &id) const
uint32_t tidGlued(const DetId &id) const