CMS 3D CMS Logo

RecHitTools.cc
Go to the documentation of this file.
2 
12 
14 
17 
18 using namespace hgcal;
19 
20 namespace {
21  template <typename DDD>
22  inline void check_ddd(const DDD* ddd) {
23  if (nullptr == ddd) {
24  throw cms::Exception("hgcal::RecHitTools") << "DDDConstants not accessibl to hgcal::RecHitTools!";
25  }
26  }
27 
28  template <typename GEOM>
29  inline void check_geom(const GEOM* geom) {
30  if (nullptr == geom) {
31  throw cms::Exception("hgcal::RecHitTools") << "Geometry not provided yet to hgcal::RecHitTools!";
32  }
33  }
34 
35  inline const HGCalDDDConstants* get_ddd(const CaloSubdetectorGeometry* geom, const HGCalDetId& detid) {
36  const HGCalGeometry* hg = static_cast<const HGCalGeometry*>(geom);
37  const HGCalDDDConstants* ddd = &(hg->topology().dddConstants());
38  check_ddd(ddd);
39  return ddd;
40  }
41 
42  inline const HGCalDDDConstants* get_ddd(const CaloSubdetectorGeometry* geom, const HGCSiliconDetId& detid) {
43  const HGCalGeometry* hg = static_cast<const HGCalGeometry*>(geom);
44  const HGCalDDDConstants* ddd = &(hg->topology().dddConstants());
45  check_ddd(ddd);
46  return ddd;
47  }
48 
49  inline const HGCalDDDConstants* get_ddd(const CaloSubdetectorGeometry* geom, const HFNoseDetId& detid) {
50  const HGCalGeometry* hg = static_cast<const HGCalGeometry*>(geom);
51  const HGCalDDDConstants* ddd = &(hg->topology().dddConstants());
52  check_ddd(ddd);
53  return ddd;
54  }
55 
56  inline const HGCalDDDConstants* get_ddd(const CaloGeometry* geom, int type, int maxLayerEE, int layer) {
57  DetId::Detector det = ((type == 0) ? DetId::Forward : ((layer > maxLayerEE) ? DetId::HGCalHSi : DetId::HGCalEE));
58  int subdet = ((type != 0) ? ForwardSubdetector::ForwardEmpty
59  : ((layer > maxLayerEE) ? ForwardSubdetector::HGCEE : ForwardSubdetector::HGCHEF));
60  const HGCalGeometry* hg = static_cast<const HGCalGeometry*>(geom->getSubdetectorGeometry(det, subdet));
61  const HGCalDDDConstants* ddd = &(hg->topology().dddConstants());
62  check_ddd(ddd);
63  return ddd;
64  }
65 
66  enum CellType {
67  CE_E_120 = 0,
68  CE_E_200 = 1,
69  CE_E_300 = 2,
70  CE_H_120 = 3,
71  CE_H_200 = 4,
72  CE_H_300 = 5,
73  CE_H_SCINT = 6,
74  EnumSize = 7
75  };
76 
77 } // namespace
78 
80  geom_ = &geom;
81  unsigned int wmaxEE(0), wmaxFH(0);
82  auto geomEE = static_cast<const HGCalGeometry*>(
84  //check if it's the new geometry
85  if (geomEE) {
86  geometryType_ = 1;
87  eeOffset_ = (geomEE->topology().dddConstants()).getLayerOffset();
88  wmaxEE = (geomEE->topology().dddConstants()).waferCount(0);
89  auto geomFH = static_cast<const HGCalGeometry*>(
91  fhOffset_ = (geomFH->topology().dddConstants()).getLayerOffset();
92  wmaxFH = (geomFH->topology().dddConstants()).waferCount(0);
93  fhLastLayer_ = fhOffset_ + (geomFH->topology().dddConstants()).lastLayer(true);
94  auto geomBH = static_cast<const HGCalGeometry*>(
96  bhOffset_ = (geomBH->topology().dddConstants()).getLayerOffset();
97  bhFirstLayer_ = bhOffset_ + (geomBH->topology().dddConstants()).firstLayer();
98  bhLastLayer_ = bhOffset_ + (geomBH->topology().dddConstants()).lastLayer(true);
99  bhMaxIphi_ = geomBH->topology().dddConstants().maxCells(true);
100  } else {
101  geometryType_ = 0;
102  geomEE =
104  eeOffset_ = (geomEE->topology().dddConstants()).getLayerOffset();
105  wmaxEE = 1 + (geomEE->topology().dddConstants()).waferMax();
106  auto geomFH =
108  fhOffset_ = (geomFH->topology().dddConstants()).getLayerOffset();
109  fhLastLayer_ = fhOffset_ + (geomFH->topology().dddConstants()).layers(true);
111  bhFirstLayer_ = bhOffset_ + 1;
112  wmaxFH = 1 + (geomFH->topology().dddConstants()).waferMax();
113  auto geomBH =
115  bhLastLayer_ = bhOffset_ + (geomBH->topology().dddConstants())->getMaxDepth(1);
116  }
117  maxNumberOfWafersPerLayer_ = std::max(wmaxEE, wmaxFH);
118  // For nose geometry
119  auto geomNose =
121  if (geomNose) {
122  maxNumberOfWafersNose_ = (geomNose->topology().dddConstants()).waferCount(0);
123  noseLastLayer_ = (geomNose->topology().dddConstants()).layers(true);
124  } else {
126  noseLastLayer_ = 0;
127  }
128 }
129 
131  DetId::Detector det = id.det();
132  int subdet = (det == DetId::HGCalEE || det == DetId::HGCalHSi || det == DetId::HGCalHSc)
134  : id.subdetId();
135  auto geom = geom_->getSubdetectorGeometry(det, subdet);
136  check_geom(geom);
137  return geom;
138 }
139 
141  auto geom = getSubdetectorGeometry(id);
143  if (id.det() == DetId::Hcal) {
144  position = geom->getGeometry(id)->getPosition();
145  } else {
146  auto hg = static_cast<const HGCalGeometry*>(geom);
147  position = hg->getPosition(id);
148  }
149  return position;
150 }
151 
152 GlobalPoint RecHitTools::getPositionLayer(int layer, bool nose) const {
153  unsigned int lay = std::abs(layer);
154  double z(0);
155  if (nose) {
156  auto geomNose =
158  if (geomNose) {
159  const HGCalDDDConstants* ddd = &(geomNose->topology().dddConstants());
160  if (ddd)
161  z = (layer > 0) ? ddd->waferZ(lay, true) : -ddd->waferZ(lay, true);
162  }
163  } else {
164  const HGCalDDDConstants* ddd = get_ddd(geom_, geometryType_, fhOffset_, lay);
165  if (geometryType_ == 1) {
166  if (lay > fhOffset_)
167  lay -= fhOffset_;
168  }
169  z = (layer > 0) ? ddd->waferZ(lay, true) : -ddd->waferZ(lay, true);
170  }
171  return GlobalPoint(0, 0, z);
172 }
173 
174 int RecHitTools::zside(const DetId& id) const {
175  int zside = 0;
176  if (id.det() == DetId::HGCalEE || id.det() == DetId::HGCalHSi) {
177  zside = HGCSiliconDetId(id).zside();
178  } else if (id.det() == DetId::HGCalHSc) {
180  } else if (id.det() == DetId::Forward && id.subdetId() == static_cast<int>(HFNose)) {
181  zside = HFNoseDetId(id).zside();
182  } else if (id.det() == DetId::Forward) {
183  zside = HGCalDetId(id).zside();
184  } else if (id.det() == DetId::Hcal && id.subdetId() == HcalEndcap) {
185  zside = HcalDetId(id).zside();
186  }
187  return zside;
188 }
189 
190 std::float_t RecHitTools::getSiThickness(const DetId& id) const {
191  auto geom = getSubdetectorGeometry(id);
192  std::float_t thick(0.37);
193  if (id.det() == DetId::HGCalEE || id.det() == DetId::HGCalHSi) {
194  const HGCSiliconDetId hid(id);
195  auto ddd = get_ddd(geom, hid);
196  thick = ddd->cellThickness(hid.layer(), hid.waferU(), hid.waferV());
197  } else if (id.det() == DetId::Forward && id.subdetId() == static_cast<int>(HFNose)) {
198  const HFNoseDetId hid(id);
199  auto ddd = get_ddd(geom, hid);
200  thick = ddd->cellThickness(hid.layer(), hid.waferU(), hid.waferV());
201  } else if (id.det() == DetId::Forward) {
202  const HGCalDetId hid(id);
203  auto ddd = get_ddd(geom, hid);
204  thick = ddd->cellThickness(hid.layer(), hid.wafer(), 0);
205  } else {
206  LogDebug("getSiThickness::InvalidSiliconDetid")
207  << "det id: " << std::hex << id.rawId() << std::dec << ":" << id.det() << " is not HGCal silicon!";
208  // It does not make any sense to return 0.37 as thickness for a detector
209  // that is not Silicon(does it?). Mark it as 0. to be able to distinguish
210  // the case.
211  thick = 0.f;
212  }
213  return thick;
214 }
215 
216 int RecHitTools::getSiThickIndex(const DetId& id) const {
217  int thickIndex = -1;
218  if (id.det() == DetId::HGCalEE || id.det() == DetId::HGCalHSi) {
219  thickIndex = HGCSiliconDetId(id).type();
220  } else if (id.det() == DetId::Forward && id.subdetId() == static_cast<int>(HFNose)) {
221  thickIndex = HFNoseDetId(id).type();
222  } else if (id.det() == DetId::Forward) {
223  float thickness = getSiThickness(id);
224  if (thickness > 99. && thickness < 121.)
225  thickIndex = 0;
226  else if (thickness > 199. && thickness < 201.)
227  thickIndex = 1;
228  else if (thickness > 299. && thickness < 301.)
229  thickIndex = 2;
230  else
231  assert(thickIndex > 0 && "ERROR - silicon thickness has a nonsensical value");
232  }
233  return thickIndex;
234 }
235 
236 std::pair<float, float> RecHitTools::getScintDEtaDPhi(const DetId& id) const {
237  if (!isScintillator(id)) {
238  LogDebug("getScintDEtaDPhi::InvalidScintDetid")
239  << "det id: " << std::hex << id.rawId() << std::dec << ":" << id.det() << " is not HGCal scintillator!";
240  return {0.f, 0.f};
241  }
242  auto cellGeom = getSubdetectorGeometry(id)->getGeometry(id);
243  return {cellGeom->etaSpan(), cellGeom->phiSpan()};
244 }
245 
246 std::float_t RecHitTools::getRadiusToSide(const DetId& id) const {
247  auto geom = getSubdetectorGeometry(id);
248  std::float_t size(std::numeric_limits<std::float_t>::max());
249  if (id.det() == DetId::HGCalEE || id.det() == DetId::HGCalHSi) {
250  const HGCSiliconDetId hid(id);
251  auto ddd = get_ddd(geom, hid);
252  size = ddd->cellSizeHex(hid.type());
253  } else if (id.det() == DetId::Forward && id.subdetId() == static_cast<int>(HFNose)) {
254  const HFNoseDetId hid(id);
255  auto ddd = get_ddd(geom, hid);
256  size = ddd->cellSizeHex(hid.type());
257  } else if (id.det() == DetId::Forward) {
258  const HGCalDetId hid(id);
259  auto ddd = get_ddd(geom, hid);
260  size = ddd->cellSizeHex(hid.waferType());
261  } else {
262  edm::LogError("getRadiusToSide::InvalidSiliconDetid")
263  << "det id: " << std::hex << id.rawId() << std::dec << ":" << id.det() << " is not HGCal silicon!";
264  }
265  return size;
266 }
267 
268 unsigned int RecHitTools::getLayer(const ForwardSubdetector type) const {
269  int layer(0);
270  switch (type) {
271  case (ForwardSubdetector::HGCEE): {
272  auto geomEE =
274  layer = (geomEE->topology().dddConstants()).layers(true);
275  break;
276  }
278  auto geomFH =
280  layer = (geomFH->topology().dddConstants()).layers(true);
281  break;
282  }
284  auto geomBH =
286  layer = (geomBH->topology().dddConstants())->getMaxDepth(1);
287  break;
288  }
290  auto geomEE =
292  layer = (geomEE->topology().dddConstants()).layers(true);
293  auto geomFH =
295  layer += (geomFH->topology().dddConstants()).layers(true);
296  auto geomBH =
298  if (geomBH)
299  layer += (geomBH->topology().dddConstants())->getMaxDepth(1);
300  auto geomBH2 =
302  if (geomBH2)
303  layer += (geomBH2->topology().dddConstants()).layers(true);
304  break;
305  }
307  auto geomHFN =
309  layer = (geomHFN->topology().dddConstants()).layers(true);
310  break;
311  }
312  default:
313  layer = 0;
314  }
315  return (unsigned int)(layer);
316 }
317 
318 unsigned int RecHitTools::getLayer(const DetId::Detector type, bool nose) const {
319  int layer;
320  switch (type) {
321  case (DetId::HGCalEE): {
322  auto geomEE =
324  layer = (geomEE->topology().dddConstants()).layers(true);
325  break;
326  }
327  case (DetId::HGCalHSi): {
328  auto geomFH =
330  layer = (geomFH->topology().dddConstants()).layers(true);
331  break;
332  }
333  case (DetId::HGCalHSc): {
334  auto geomBH =
336  layer = (geomBH->topology().dddConstants()).layers(true);
337  break;
338  }
339  case (DetId::Forward): {
340  if (nose) {
341  auto geomHFN = static_cast<const HGCalGeometry*>(
343  layer = (geomHFN->topology().dddConstants()).layers(true);
344  } else {
345  auto geomEE = static_cast<const HGCalGeometry*>(
347  layer = (geomEE->topology().dddConstants()).layers(true);
348  auto geomFH = static_cast<const HGCalGeometry*>(
350  layer += (geomFH->topology().dddConstants()).layers(true);
351  }
352  break;
353  }
354  default:
355  layer = 0;
356  }
357  return (unsigned int)(layer);
358 }
359 
360 unsigned int RecHitTools::getLayer(const DetId& id) const {
361  unsigned int layer = std::numeric_limits<unsigned int>::max();
362  if (id.det() == DetId::HGCalEE || id.det() == DetId::HGCalHSi) {
363  layer = HGCSiliconDetId(id).layer();
364  } else if (id.det() == DetId::HGCalHSc) {
365  layer = HGCScintillatorDetId(id).layer();
366  } else if (id.det() == DetId::Forward && id.subdetId() == static_cast<int>(HFNose)) {
367  layer = HFNoseDetId(id).layer();
368  } else if (id.det() == DetId::Forward) {
369  layer = HGCalDetId(id).layer();
370  } else if (id.det() == DetId::Hcal && id.subdetId() == HcalEndcap) {
371  layer = HcalDetId(id).depth();
372  }
373  return layer;
374 }
375 
376 unsigned int RecHitTools::getLayerWithOffset(const DetId& id) const {
377  unsigned int layer = getLayer(id);
378  if (id.det() == DetId::Forward && id.subdetId() == HGCHEF) {
379  layer += fhOffset_;
380  } else if (id.det() == DetId::HGCalHSi || id.det() == DetId::HGCalHSc) {
381  // DetId::HGCalHSc hits include the offset w.r.t. EE already
382  layer += fhOffset_;
383  } else if (id.det() == DetId::Hcal && id.subdetId() == HcalEndcap) {
384  layer += bhOffset_;
385  }
386  // no need to add offset for HFnose
387  return layer;
388 }
389 
390 std::pair<int, int> RecHitTools::getWafer(const DetId& id) const {
392  int waferV = 0;
393  if ((id.det() == DetId::HGCalEE) || (id.det() == DetId::HGCalHSi)) {
394  waferU = HGCSiliconDetId(id).waferU();
395  waferV = HGCSiliconDetId(id).waferV();
396  } else if (id.det() == DetId::Forward && id.subdetId() == static_cast<int>(HFNose)) {
397  waferU = HFNoseDetId(id).waferU();
398  waferV = HFNoseDetId(id).waferV();
399  } else if (id.det() == DetId::Forward) {
400  waferU = HGCalDetId(id).wafer();
401  } else {
402  edm::LogError("getWafer::InvalidSiliconDetid")
403  << "det id: " << std::hex << id.rawId() << std::dec << ":" << id.det() << " is not HGCal silicon!";
404  }
405  return std::pair<int, int>(waferU, waferV);
406 }
407 
408 std::pair<int, int> RecHitTools::getCell(const DetId& id) const {
409  int cellU = std::numeric_limits<int>::max();
410  int cellV = 0;
411  if ((id.det() == DetId::HGCalEE) || (id.det() == DetId::HGCalHSi)) {
412  cellU = HGCSiliconDetId(id).cellU();
413  cellV = HGCSiliconDetId(id).cellV();
414  } else if (id.det() == DetId::Forward && id.subdetId() == static_cast<int>(HFNose)) {
415  cellU = HFNoseDetId(id).cellU();
416  cellV = HFNoseDetId(id).cellV();
417  } else if (id.det() == DetId::Forward) {
418  cellU = HGCalDetId(id).cell();
419  } else {
420  edm::LogError("getCell::InvalidSiliconDetid")
421  << "det id: " << std::hex << id.rawId() << std::dec << ":" << id.det() << " is not HGCal silicon!";
422  }
423  return std::pair<int, int>(cellU, cellV);
424 }
425 
426 bool RecHitTools::isHalfCell(const DetId& id) const {
427  bool ishalf = false;
428  if (id.det() == DetId::Forward) {
429  HGCalDetId hid(id);
430  auto geom = getSubdetectorGeometry(hid);
431  auto ddd = get_ddd(geom, hid);
432  const int waferType = ddd->waferTypeT(hid.waferType());
433  return ddd->isHalfCell(waferType, hid.cell());
434  }
435  //new geometry is always false
436  return ishalf;
437 }
438 
439 int RecHitTools::getCellType(const DetId& id) const {
440  auto layer_number = getLayerWithOffset(id);
441  auto thickness = getSiThickIndex(id);
442  auto geomNose =
444  auto isNose = geomNose ? true : false;
445  auto isEELayer = (layer_number <= lastLayerEE(isNose));
446  auto isScint = isScintillator(id);
447  int layerType = -1;
448 
449  if (isScint) {
450  layerType = CE_H_SCINT;
451  }
452  if (isEELayer) {
453  if (thickness == 0) {
454  layerType = CE_E_120;
455  } else if (thickness == 1) {
456  layerType = CE_E_200;
457  } else if (thickness == 2) {
458  layerType = CE_E_300;
459  }
460  } else {
461  if (thickness == 0) {
462  layerType = CE_H_120;
463  } else if (thickness == 1) {
464  layerType = CE_H_200;
465  } else if (thickness == 2) {
466  layerType = CE_H_300;
467  }
468  }
469  assert(layerType != -1);
470  return layerType;
471 }
472 
473 bool RecHitTools::isSilicon(const DetId& id) const {
474  return (id.det() == DetId::HGCalEE || id.det() == DetId::HGCalHSi ||
475  (id.det() == DetId::Forward && id.subdetId() == static_cast<int>(HFNose)));
476 }
477 
478 bool RecHitTools::isScintillator(const DetId& id) const { return id.det() == DetId::HGCalHSc; }
479 
480 bool RecHitTools::isOnlySilicon(const unsigned int layer) const {
481  // HFnose TODO
482  bool isonlysilicon = (layer % bhLastLayer_) < bhOffset_;
483  return isonlysilicon;
484 }
485 
486 float RecHitTools::getEta(const GlobalPoint& position, const float& vertex_z) const {
487  GlobalPoint corrected_position = GlobalPoint(position.x(), position.y(), position.z() - vertex_z);
488  return corrected_position.eta();
489 }
490 
491 float RecHitTools::getEta(const DetId& id, const float& vertex_z) const {
493  float eta = getEta(position, vertex_z);
494  return eta;
495 }
496 
498  float phi = atan2(position.y(), position.x());
499  return phi;
500 }
501 
502 float RecHitTools::getPhi(const DetId& id) const {
504  float phi = atan2(position.y(), position.x());
505  return phi;
506 }
507 
508 float RecHitTools::getPt(const GlobalPoint& position, const float& hitEnergy, const float& vertex_z) const {
509  float eta = getEta(position, vertex_z);
510  float pt = hitEnergy / cosh(eta);
511  return pt;
512 }
513 
514 float RecHitTools::getPt(const DetId& id, const float& hitEnergy, const float& vertex_z) const {
516  float eta = getEta(position, vertex_z);
517  float pt = hitEnergy / cosh(eta);
518  return pt;
519 }
520 
521 std::pair<uint32_t, uint32_t> RecHitTools::firstAndLastLayer(DetId::Detector det, int subdet) const {
522  if ((det == DetId::HGCalEE) || ((det == DetId::Forward) && (subdet == HGCEE))) {
523  return std::make_pair(eeOffset_ + 1, fhOffset_);
524  } else if ((det == DetId::HGCalHSi) || ((det == DetId::Forward) && (subdet == HGCHEF))) {
525  return std::make_pair(fhOffset_ + 1, fhLastLayer_);
526  } else if ((det == DetId::Forward) && (subdet == HFNose)) {
527  return std::make_pair(1, noseLastLayer_);
528  } else {
529  return std::make_pair(bhFirstLayer_, bhLastLayer_);
530  }
531 }
532 
533 bool RecHitTools::maskCell(const DetId& id, int corners) const {
534  if (id.det() == DetId::Hcal) {
535  return false;
536  } else {
537  auto hg = static_cast<const HGCalGeometry*>(getSubdetectorGeometry(id));
538  return hg->topology().dddConstants().maskCell(id, corners);
539  }
540 }
bool maskCell(const DetId &id, int corners) const
bool isHalfCell(const DetId &) const
Definition: RecHitTools.cc:426
double waferZ(int layer, bool reco) const
std::pair< int, int > getWafer(const DetId &) const
Definition: RecHitTools.cc:390
unsigned int bhFirstLayer_
Definition: RecHitTools.h:91
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
Definition: RecHitTools.cc:130
int wafer() const
get the wafer #
Definition: HGCalDetId.h:40
unsigned int fhLastLayer_
Definition: RecHitTools.h:91
int waferType() const
get the wafer type
Definition: HGCalDetId.h:43
int layer() const
get the layer #
Definition: HFNoseDetId.h:57
unsigned int getLayer(DetId::Detector type, bool nose=false) const
Definition: RecHitTools.cc:318
constexpr int zside() const
get the z-side of the cell (1/-1)
Definition: HcalDetId.h:141
int32_t waferU(const int32_t index)
float getPhi(const GlobalPoint &position) const
Definition: RecHitTools.cc:497
unsigned int bhLastLayer_
Definition: RecHitTools.h:91
int waferU() const
Definition: HFNoseDetId.h:76
bool maskCell(const DetId &id, int corners=3) const
Definition: RecHitTools.cc:533
T eta() const
Definition: PV3DBase.h:73
unsigned int maxNumberOfWafersNose_
Definition: RecHitTools.h:92
int getCellType(const DetId &id) const
Definition: RecHitTools.cc:439
std::float_t getRadiusToSide(const DetId &) const
Definition: RecHitTools.cc:246
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
int waferTypeT(int wafer) const
unsigned int fhOffset_
Definition: RecHitTools.h:91
constexpr int zside() const
get the z-side of the cell (1/-1)
Log< level::Error, false > LogError
CellType
Definition: RecHitTools.cc:66
std::pair< float, float > getScintDEtaDPhi(const DetId &) const
Definition: RecHitTools.cc:236
unsigned int bhOffset_
Definition: RecHitTools.h:91
ForwardSubdetector
assert(be >=bs)
int zside(const DetId &id) const
Definition: RecHitTools.cc:174
constexpr int32_t waferV() const
bool isHalfCell(int waferType, int cell) const
unsigned int maxNumberOfWafersPerLayer_
Definition: RecHitTools.h:92
GlobalPoint getPosition(const DetId &id) const
Definition: RecHitTools.cc:140
constexpr int32_t zside() const
get the z-side of the cell (1/-1)
const CaloGeometry * geom_
Definition: RecHitTools.h:90
bool isSilicon(const DetId &) const
Definition: RecHitTools.cc:473
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int zside() const
get the z-side of the cell (1/-1)
Definition: HGCalDetId.h:49
int cellU() const
get the cell #&#39;s in u,v or in x,y
Definition: HFNoseDetId.h:60
bool isOnlySilicon(const unsigned int layer) const
Definition: RecHitTools.cc:480
constexpr int32_t cellU() const
get the cell #&#39;s in u,v or in x,y
int type() const
get the type
Definition: HFNoseDetId.h:51
const HGCalTopology & topology() const
int cell() const
get the absolute value of the cell #&#39;s in x and y
Definition: HGCalDetId.h:37
double cellSizeHex(int type) const
virtual std::shared_ptr< const CaloCellGeometry > getGeometry(const DetId &id) const
Get the cell geometry of a given detector id. Should return false if not found.
std::pair< uint32_t, uint32_t > firstAndLastLayer(DetId::Detector det, int subdet) const
Definition: RecHitTools.cc:521
Definition: DetId.h:17
unsigned int eeOffset_
Definition: RecHitTools.h:91
int cellV() const
Definition: HFNoseDetId.h:61
Detector
Definition: DetId.h:24
std::pair< int, int > getCell(const DetId &) const
Definition: RecHitTools.cc:408
double cellThickness(int layer, int waferU, int waferV) const
int layer() const
get the layer #
Definition: HGCalDetId.h:46
float getEta(const GlobalPoint &position, const float &vertex_z=0.) const
Definition: RecHitTools.cc:486
GlobalPoint getPositionLayer(int layer, bool nose=false) const
Definition: RecHitTools.cc:152
int waferV() const
Definition: HFNoseDetId.h:79
constexpr int32_t layer() const
get the layer #
constexpr int32_t cellV() const
unsigned int noseLastLayer_
Definition: RecHitTools.h:91
void setGeometry(CaloGeometry const &)
Definition: RecHitTools.cc:79
bool isScintillator(const DetId &) const
Definition: RecHitTools.cc:478
GlobalPoint getPosition(const DetId &id, bool debug=false) const
constexpr int layer() const
get the layer #
static int position[264][3]
Definition: ReadPGInfo.cc:289
int32_t waferV(const int32_t index)
constexpr int32_t waferU() const
constexpr int32_t type() const
get the type
int getSiThickIndex(const DetId &) const
Definition: RecHitTools.cc:216
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:34
float getPt(const GlobalPoint &position, const float &hitEnergy, const float &vertex_z=0.) const
Definition: RecHitTools.cc:508
const HGCalDDDConstants & dddConstants() const
Definition: HGCalTopology.h:98
unsigned int lastLayerEE(bool nose=false) const
Definition: RecHitTools.h:76
std::float_t getSiThickness(const DetId &) const
Definition: RecHitTools.cc:190
int zside() const
get the z-side of the cell (1/-1)
Definition: HFNoseDetId.h:54
#define LogDebug(id)
unsigned int lastLayer(bool nose=false) const
Definition: RecHitTools.h:80
unsigned int getLayerWithOffset(const DetId &) const
Definition: RecHitTools.cc:376
constexpr int depth() const
get the tower depth
Definition: HcalDetId.h:164