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")
25  << "DDDConstants not accessibl to hgcal::RecHitTools!";
26  }
27  }
28 
29  template<typename GEOM>
30  inline void check_geom(const GEOM* geom) {
31  if( nullptr == geom ) {
32  throw cms::Exception("hgcal::RecHitTools")
33  << "Geometry not provided yet to hgcal::RecHitTools!";
34  }
35  }
36 
37  inline const HGCalDDDConstants* get_ddd(const CaloSubdetectorGeometry* geom,
38  const HGCalDetId& detid) {
39  const HGCalGeometry* hg = static_cast<const HGCalGeometry*>(geom);
40  const HGCalDDDConstants* ddd = &(hg->topology().dddConstants());
41  check_ddd(ddd);
42  return ddd;
43  }
44 
45  inline const HGCalDDDConstants* get_ddd(const CaloSubdetectorGeometry* geom,
46  const HGCSiliconDetId& detid) {
47  const HGCalGeometry* hg = static_cast<const HGCalGeometry*>(geom);
48  const HGCalDDDConstants* ddd = &(hg->topology().dddConstants());
49  check_ddd(ddd);
50  return ddd;
51  }
52 
53  inline const HGCalDDDConstants* get_ddd(const CaloGeometry* geom,
54  int type, int maxLayerEE,
55  int layer) {
56  DetId::Detector det = ((type == 0) ? DetId::Forward :
57  ((layer > maxLayerEE) ? DetId::HGCalHSi :
59  int subdet = ((type != 0) ? ForwardSubdetector::ForwardEmpty :
60  ((layer > maxLayerEE) ? ForwardSubdetector::HGCEE :
62  const HGCalGeometry* hg = static_cast<const HGCalGeometry*>(geom->getSubdetectorGeometry(det,subdet));
63  const HGCalDDDConstants* ddd = &(hg->topology().dddConstants());
64  check_ddd(ddd);
65  return ddd;
66  }
67 
68 }
69 
71 }
72 
74 
76  es.get<CaloGeometryRecord>().get(geom);
77 
78  geom_ = geom.product();
79  unsigned int wmaxEE(0), wmaxFH(0);
81  //check if it's the new geometry
82  if(geomEE) {
83  geometryType_ = 1;
84  fhOffset_ = (geomEE->topology().dddConstants()).layers(true);
85  wmaxEE = (geomEE->topology().dddConstants()).waferCount(0);
88  wmaxFH = (geomFH->topology().dddConstants()).waferCount(0);
89  fhLastLayer_ = fhOffset_ + (geomFH->topology().dddConstants()).layers(true);
90  }
91  else {
92  geometryType_ = 0;
94  fhOffset_ = (geomEE->topology().dddConstants()).layers(true);
95  wmaxEE = 1 + (geomEE->topology().dddConstants()).waferMax();
97  bhOffset_ = fhOffset_ + (geomFH->topology().dddConstants()).layers(true);
98  wmaxFH = 1 + (geomFH->topology().dddConstants()).waferMax();
100  }
101  maxNumberOfWafersPerLayer_ = std::max(wmaxEE,wmaxFH);
102 }
103 
105  DetId::Detector det = id.det();
106  int subdet = (det == DetId::HGCalEE || det == DetId::HGCalHSi || det == DetId::HGCalHSc) ? ForwardSubdetector::ForwardEmpty : id.subdetId();
107  auto geom = geom_->getSubdetectorGeometry(det,subdet);
108  check_geom(geom);
109  return geom;
110 }
111 
113  auto geom = getSubdetectorGeometry(id);
115  if (id.det() == DetId::Hcal) {
116  position = geom->getGeometry(id)->getPosition();
117  } else {
118  auto hg = static_cast<const HGCalGeometry*>(geom);
119  position = hg->getPosition(id);
120  }
121  return position;
122 }
123 
125  int lay = std::abs(layer);
126  const HGCalDDDConstants* ddd = get_ddd(geom_, geometryType_, fhOffset_, lay);
127  double z = (layer > 0) ? ddd->waferZ(lay,true) : -ddd->waferZ(lay,true);
128  return GlobalPoint(0,0,z);
129 }
130 
131 int RecHitTools::zside(const DetId& id) const {
132  int zside = 0;
133  if (id.det() == DetId::Forward) {
134  zside = HGCalDetId(id).zside();
135  } else if( id.det() == DetId::Hcal && id.subdetId() == HcalEndcap) {
136  zside = HcalDetId(id).zside();
137  } else if (id.det() == DetId::HGCalEE || id.det() == DetId::HGCalHSi) {
138  zside = HGCSiliconDetId(id).zside();
139  } else if (id.det() == DetId::HGCalHSc) {
140  zside = HGCScintillatorDetId(id).zside();
141  }
142  return zside;
143 }
144 
145 std::float_t RecHitTools::getSiThickness(const DetId& id) const {
146  auto geom = getSubdetectorGeometry(id);
147  std::float_t thick(0.37);
148  if (id.det() == DetId::Forward) {
149  const HGCalDetId hid(id);
150  auto ddd = get_ddd(geom,hid);
151  thick = ddd->cellThickness(hid.layer(),hid.wafer(),0);
152  } else if (id.det() == DetId::HGCalEE || id.det() == DetId::HGCalHSi) {
153  const HGCSiliconDetId hid(id);
154  auto ddd = get_ddd(geom,hid);
155  thick = ddd->cellThickness(hid.layer(),hid.waferU(),hid.waferV());
156  } else {
157  LogDebug("getSiThickness::InvalidSiliconDetid")
158  << "det id: " << std::hex << id.rawId() << std::dec << ":"
159  << id.det() << " is not HGCal silicon!";
160  // It does not make any sense to return 0.37 as thickness for a detector
161  // that is not Silicon(does it?). Mark it as 0. to be able to distinguish
162  // the case.
163  thick = 0.f;
164  }
165  return thick;
166 }
167 
168 int RecHitTools::getSiThickIndex(const DetId& id) const {
169  int thickIndex = -1;
170  if (id.det() == DetId::Forward) {
171  float thickness = getSiThickness(id);
172  if (thickness > 99. && thickness < 121.)
173  thickIndex = 0;
174  else if (thickness > 199. && thickness < 201.)
175  thickIndex = 1;
176  else if (thickness > 299. && thickness < 301.)
177  thickIndex = 2;
178  else
179  assert(thickIndex > 0 && "ERROR - silicon thickness has a nonsensical value");
180  } else if (id.det() == DetId::HGCalEE || id.det() == DetId::HGCalHSi) {
181  thickIndex = HGCSiliconDetId(id).type();
182  }
183  return thickIndex;
184 }
185 
186 std::float_t RecHitTools::getRadiusToSide(const DetId& id) const {
187  auto geom = getSubdetectorGeometry(id);
189  if (id.det() == DetId::Forward) {
190  const HGCalDetId hid(id);
191  auto ddd = get_ddd(geom,hid);
192  size = ddd->cellSizeHex(hid.waferType());
193  } else if (id.det() == DetId::HGCalEE || id.det() == DetId::HGCalHSi) {
194  const HGCSiliconDetId hid(id);
195  auto ddd = get_ddd(geom,hid);
196  size = ddd->cellSizeHex(hid.type());
197  } else {
198  edm::LogError("getRadiusToSide::InvalidSiliconDetid")
199  << "det id: " << std::hex << id.rawId() << std::dec << ":"
200  << id.det() << " is not HGCal silicon!";
201  }
202  return size;
203 }
204 
205 unsigned int RecHitTools::getLayer(const ForwardSubdetector type) const {
206 
207  int layer;
208  switch (type) {
211  layer = (geomEE->topology().dddConstants()).layers(true);
212  break;
213  }
216  layer = (geomFH->topology().dddConstants()).layers(true);
217  break;
218  }
221  layer = (geomBH->topology().dddConstants())->getMaxDepth(1);
222  break;
223  }
226  layer = (geomEE->topology().dddConstants()).layers(true);
228  layer += (geomFH->topology().dddConstants()).layers(true);
230  layer += (geomBH->topology().dddConstants())->getMaxDepth(1);
231  break;
232  }
233  default: layer = 0;
234  }
235  return (unsigned int)(layer);
236 }
237 
238 unsigned int RecHitTools::getLayer(const DetId::Detector type) const {
239 
240  int layer;
241  switch (type) {
242  case(DetId::HGCalEE): {
243  auto geomEE = static_cast<const HGCalGeometry*>(geom_->getSubdetectorGeometry(type,ForwardSubdetector::ForwardEmpty));
244  layer = (geomEE->topology().dddConstants()).layers(true);
245  break;
246  }
247  case (DetId::HGCalHSi): {
248  auto geomFH = static_cast<const HGCalGeometry*>(geom_->getSubdetectorGeometry(type,ForwardSubdetector::ForwardEmpty));
249  layer = (geomFH->topology().dddConstants()).layers(true);
250  break;
251  }
252  case (DetId::HGCalHSc): {
253  auto geomBH = static_cast<const HGCalGeometry*>(geom_->getSubdetectorGeometry(type,ForwardSubdetector::ForwardEmpty));
254  layer = (geomBH->topology().dddConstants()).layers(true);
255  break;
256  }
257  case (DetId::Forward): {
259  layer = (geomEE->topology().dddConstants()).layers(true);
261  layer += (geomFH->topology().dddConstants()).layers(true);
262  break;
263  }
264  default: layer = 0;
265  }
266  return (unsigned int)(layer);
267 }
268 
269 unsigned int RecHitTools::getLayer(const DetId& id) const {
270  unsigned int layer = std::numeric_limits<unsigned int>::max();
271  if (id.det() == DetId::Forward) {
272  layer = HGCalDetId(id).layer();
273  } else if (id.det() == DetId::Hcal && id.subdetId() == HcalEndcap) {
274  layer = HcalDetId(id).depth();
275  } else if (id.det() == DetId::HGCalEE || id.det() == DetId::HGCalHSi) {
276  layer = HGCSiliconDetId(id).layer();
277  } else if (id.det() == DetId::HGCalHSc) {
278  layer = HGCScintillatorDetId(id).layer();
279  }
280  return layer;
281 }
282 
283 unsigned int RecHitTools::getLayerWithOffset(const DetId& id) const {
284  unsigned int layer = getLayer(id);
285  if (id.det() == DetId::Forward && id.subdetId() == HGCHEF ) {
286  layer += fhOffset_;
287  } else if (id.det() == DetId::HGCalHSi || id.det() == DetId::HGCalHSc) {
288  layer += fhOffset_;
289  } else if (id.det() == DetId::Hcal && id.subdetId() == HcalEndcap) {
290  layer += bhOffset_;
291  }
292  return layer;
293 }
294 
295 std::pair<int,int> RecHitTools::getWafer(const DetId& id) const {
296  int waferU = std::numeric_limits<int>::max();
297  int waferV = 0;
298  if (id.det() == DetId::Forward) {
299  waferU = HGCalDetId(id).wafer();
300  } else if ((id.det() == DetId::HGCalEE) || (id.det() == DetId::HGCalHSi)) {
301  waferU = HGCSiliconDetId(id).waferU();
302  waferV = HGCSiliconDetId(id).waferV();
303  } else {
304  edm::LogError("getWafer::InvalidSiliconDetid")
305  << "det id: " << std::hex << id.rawId() << std::dec << ":"
306  << id.det() << " is not HGCal silicon!";
307  }
308  return std::pair<int,int>(waferU,waferV);
309 }
310 
311 std::pair<int,int> RecHitTools::getCell(const DetId& id) const {
312  int cellU = std::numeric_limits<int>::max();
313  int cellV = 0;
314  if (id.det() == DetId::Forward) {
315  cellU = HGCalDetId(id).cell();
316  } else if ((id.det() == DetId::HGCalEE) || (id.det() == DetId::HGCalHSi)) {
317  cellU = HGCSiliconDetId(id).cellU();
318  cellV = HGCSiliconDetId(id).cellV();
319  } else {
320  edm::LogError("getCell::InvalidSiliconDetid")
321  << "det id: " << std::hex << id.rawId() << std::dec << ":"
322  << id.det() << " is not HGCal silicon!";
323  }
324  return std::pair<int,int>(cellU,cellV);
325 }
326 
327 bool RecHitTools::isHalfCell(const DetId& id) const {
328  bool ishalf = false;
329  if (id.det() == DetId::Forward) {
330  HGCalDetId hid(id);
331  auto geom = getSubdetectorGeometry(hid);
332  auto ddd = get_ddd(geom,hid);
333  const int waferType = ddd->waferTypeT(hid.waferType());
334  return ddd->isHalfCell(waferType,hid.cell());
335  }
336  //new geometry is always false
337  return ishalf;
338 }
339 
340 float RecHitTools::getEta(const GlobalPoint& position, const float& vertex_z) const {
341  GlobalPoint corrected_position = GlobalPoint(position.x(), position.y(), position.z()-vertex_z);
342  return corrected_position.eta();
343 }
344 
345 float RecHitTools::getEta(const DetId& id, const float& vertex_z) const {
347  float eta = getEta(position, vertex_z);
348  return eta;
349 }
350 
352  float phi = atan2(position.y(),position.x());
353  return phi;
354 }
355 
356 float RecHitTools::getPhi(const DetId& id) const {
358  float phi = atan2(position.y(),position.x());
359  return phi;
360 }
361 
362 float RecHitTools::getPt(const GlobalPoint& position, const float& hitEnergy, const float& vertex_z) const {
363  float eta = getEta(position, vertex_z);
364  float pt = hitEnergy / cosh(eta);
365  return pt;
366 }
367 
368 float RecHitTools::getPt(const DetId& id, const float& hitEnergy, const float& vertex_z) const {
370  float eta = getEta(position, vertex_z);
371  float pt = hitEnergy / cosh(eta);
372  return pt;
373 }
374 
375 bool RecHitTools::maskCell(const DetId& id, int corners) const {
376  if (id.det() == DetId::Hcal) {
377  return false;
378  } else {
379  auto hg = static_cast<const HGCalGeometry*>(getSubdetectorGeometry(id));
380  return hg->topology().dddConstants().maskCell(id, corners);
381  }
382 }
#define LogDebug(id)
size
Write out results.
bool isHalfCell(int waferType, int cell) const
GlobalPoint getPositionLayer(int layer) const
Definition: RecHitTools.cc:124
type
Definition: HCALResponse.h:21
bool isHalfCell(const DetId &) const
Definition: RecHitTools.cc:327
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:49
void getEvent(const edm::Event &)
Definition: RecHitTools.cc:70
unsigned int fhLastLayer_
Definition: RecHitTools.h:65
std::vector< LayerSetAndLayers > layers(const SeedingLayerSetsHits &sets)
Definition: LayerTriplets.cc:4
int zside(const DetId &id) const
Definition: RecHitTools.cc:131
double cellSizeHex(int type) const
double cellThickness(int layer, int waferU, int waferV) const
int zside() const
get the z-side of the cell (1/-1)
Definition: HcalDetId.h:149
int waferU() const
int cellV() const
unsigned int getLayer(DetId::Detector type) const
Definition: RecHitTools.cc:238
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
T y() const
Definition: PV3DBase.h:63
unsigned int fhOffset_
Definition: RecHitTools.h:65
GlobalPoint getPosition(const DetId &id) const
bool ev
int zside() const
get the z-side of the cell (1/-1)
std::pair< int, int > getWafer(const DetId &) const
Definition: RecHitTools.cc:295
unsigned int bhOffset_
Definition: RecHitTools.h:65
ForwardSubdetector
int cellU() const
get the cell #&#39;s in u,v or in x,y
float getEta(const GlobalPoint &position, const float &vertex_z=0.) const
Definition: RecHitTools.cc:340
void getEventSetup(const edm::EventSetup &)
Definition: RecHitTools.cc:73
int depth() const
get the tower depth
Definition: HcalDetId.h:166
float getPhi(const GlobalPoint &position) const
Definition: RecHitTools.cc:351
float getPt(const GlobalPoint &position, const float &hitEnergy, const float &vertex_z=0.) const
Definition: RecHitTools.cc:362
int zside() const
get the z-side of the cell (1/-1)
Definition: HGCalDetId.h:51
unsigned int maxNumberOfWafersPerLayer_
Definition: RecHitTools.h:65
int type() const
get the type
int layer() const
get the layer #
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
Definition: RecHitTools.cc:104
T z() const
Definition: PV3DBase.h:64
const HGCalTopology & topology() const
int waferV() const
std::float_t getRadiusToSide(const DetId &) const
Definition: RecHitTools.cc:186
const CaloGeometry * geom_
Definition: RecHitTools.h:64
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int wafer() const
get the wafer #
Definition: HGCalDetId.h:42
double waferZ(int layer, bool reco) const
std::float_t getSiThickness(const DetId &) const
Definition: RecHitTools.cc:145
Definition: DetId.h:18
int cell() const
get the absolute value of the cell #&#39;s in x and y
Definition: HGCalDetId.h:39
unsigned int getLayerWithOffset(const DetId &) const
Definition: RecHitTools.cc:283
const HGCalDDDConstants & dddConstants() const
int getSiThickIndex(const DetId &) const
Definition: RecHitTools.cc:168
Detector
Definition: DetId.h:26
bool maskCell(const DetId &id, int corners) const
GlobalPoint getPosition(const DetId &id) const
Definition: RecHitTools.cc:112
int layer() const
get the layer #
T eta() const
Definition: PV3DBase.h:76
int zside() const
get the z-side of the cell (1/-1)
static int position[264][3]
Definition: ReadPGInfo.cc:509
T get() const
Definition: EventSetup.h:71
std::shared_ptr< const CaloCellGeometry > getGeometry(const DetId &id) const
Get the cell geometry of a given detector id.
Definition: CaloGeometry.cc:85
int waferType() const
get the wafer type
Definition: HGCalDetId.h:45
T x() const
Definition: PV3DBase.h:62
int waferTypeT(int wafer) const
T const * product() const
Definition: ESHandle.h:86
int layer() const
get the layer #
Definition: HGCalDetId.h:48
std::pair< int, int > getCell(const DetId &) const
Definition: RecHitTools.cc:311
bool maskCell(const DetId &id, int corners=3) const
Definition: RecHitTools.cc:375