CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
hgcal::RecHitTools Class Reference

#include <RecHitTools.h>

Public Member Functions

unsigned int firstLayerBH () const
 
std::pair< int, int > getCell (const DetId &) const
 
float getEta (const GlobalPoint &position, const float &vertex_z=0.) const
 
float getEta (const DetId &id, const float &vertex_z=0.) const
 
void getEvent (const edm::Event &)
 
void getEventSetup (const edm::EventSetup &)
 
const CaloGeometrygetGeometry () const
 
int getGeometryType () const
 
unsigned int getLayer (DetId::Detector type, bool nose=false) const
 
unsigned int getLayer (ForwardSubdetector type) const
 
unsigned int getLayer (const DetId &) const
 
unsigned int getLayerWithOffset (const DetId &) const
 
float getPhi (const GlobalPoint &position) const
 
float getPhi (const DetId &id) const
 
GlobalPoint getPosition (const DetId &id) const
 
GlobalPoint getPositionLayer (int layer, bool nose=false) const
 
float getPt (const GlobalPoint &position, const float &hitEnergy, const float &vertex_z=0.) const
 
float getPt (const DetId &id, const float &hitEnergy, const float &vertex_z=0.) const
 
std::float_t getRadiusToSide (const DetId &) const
 
int getScintMaxIphi () const
 
int getSiThickIndex (const DetId &) const
 
std::float_t getSiThickness (const DetId &) const
 
const CaloSubdetectorGeometrygetSubdetectorGeometry (const DetId &id) const
 
std::pair< int, int > getWafer (const DetId &) const
 
bool isHalfCell (const DetId &) const
 
bool isOnlySilicon (const unsigned int layer) const
 
bool isSilicon (const DetId &) const
 
unsigned int lastLayerBH () const
 
unsigned int lastLayerEE () const
 
unsigned int lastLayerFH () const
 
bool maskCell (const DetId &id, int corners=3) const
 
unsigned int maxNumberOfWafersPerLayer (bool nose=false) const
 
 RecHitTools ()
 
int zside (const DetId &id) const
 
 ~RecHitTools ()
 

Private Attributes

unsigned int bhLastLayer_
 
int bhMaxIphi_
 
unsigned int bhOffset_
 
unsigned int fhLastLayer_
 
unsigned int fhOffset_
 
const CaloGeometrygeom_
 
int geometryType_
 
unsigned int maxNumberOfWafersNose_
 
unsigned int maxNumberOfWafersPerLayer_
 

Detailed Description

Definition at line 20 of file RecHitTools.h.

Constructor & Destructor Documentation

hgcal::RecHitTools::RecHitTools ( )
inline

Definition at line 22 of file RecHitTools.h.

22 : geom_(nullptr), fhOffset_(0), bhOffset_(0), fhLastLayer_(0), geometryType_(0) {}
unsigned int fhLastLayer_
Definition: RecHitTools.h:75
unsigned int fhOffset_
Definition: RecHitTools.h:75
unsigned int bhOffset_
Definition: RecHitTools.h:75
const CaloGeometry * geom_
Definition: RecHitTools.h:74
hgcal::RecHitTools::~RecHitTools ( )
inline

Definition at line 23 of file RecHitTools.h.

References position, and ecaldqm::zside().

23 {}

Member Function Documentation

unsigned int hgcal::RecHitTools::firstLayerBH ( ) const
inline

Definition at line 64 of file RecHitTools.h.

Referenced by HGCal3DClustering::makeClusters(), and HGCalDepthPreClusterer::makePreClusters().

64 { return bhOffset_ + 1; }
unsigned int bhOffset_
Definition: RecHitTools.h:75
std::pair< int, int > RecHitTools::getCell ( const DetId id) const

Definition at line 386 of file RecHitTools.cc.

References HGCalDetId::cell(), HFNoseDetId::cellU(), HGCSiliconDetId::cellU(), HFNoseDetId::cellV(), HGCSiliconDetId::cellV(), TauDecayModes::dec, DetId::Forward, HFNose, DetId::HGCalEE, DetId::HGCalHSi, and SiStripPI::max.

386  {
387  int cellU = std::numeric_limits<int>::max();
388  int cellV = 0;
389  if ((id.det() == DetId::HGCalEE) || (id.det() == DetId::HGCalHSi)) {
390  cellU = HGCSiliconDetId(id).cellU();
391  cellV = HGCSiliconDetId(id).cellV();
392  } else if (id.det() == DetId::Forward && id.subdetId() == static_cast<int>(HFNose)) {
393  cellU = HFNoseDetId(id).cellU();
394  cellV = HFNoseDetId(id).cellV();
395  } else if (id.det() == DetId::Forward) {
396  cellU = HGCalDetId(id).cell();
397  } else {
398  edm::LogError("getCell::InvalidSiliconDetid")
399  << "det id: " << std::hex << id.rawId() << std::dec << ":" << id.det() << " is not HGCal silicon!";
400  }
401  return std::pair<int, int>(cellU, cellV);
402 }
int cellU() const
get the cell #&#39;s in u,v or in x,y
Definition: HFNoseDetId.h:58
int cellV() const
int cellU() const
get the cell #&#39;s in u,v or in x,y
int cellV() const
Definition: HFNoseDetId.h:59
int cell() const
get the absolute value of the cell #&#39;s in x and y
Definition: HGCalDetId.h:37
float RecHitTools::getEta ( const GlobalPoint position,
const float &  vertex_z = 0. 
) const

Definition at line 429 of file RecHitTools.cc.

References PV3DBase< T, PVType, FrameType >::eta(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by getEta(), getPt(), and HGCalIsoCalculator::setRecHits().

429  {
430  GlobalPoint corrected_position = GlobalPoint(position.x(), position.y(), position.z() - vertex_z);
431  return corrected_position.eta();
432 }
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
T y() const
Definition: PV3DBase.h:60
T z() const
Definition: PV3DBase.h:61
T eta() const
Definition: PV3DBase.h:73
T x() const
Definition: PV3DBase.h:59
float RecHitTools::getEta ( const DetId id,
const float &  vertex_z = 0. 
) const

Definition at line 434 of file RecHitTools.cc.

References PVValHelper::eta, getEta(), getPosition(), and position.

434  {
436  float eta = getEta(position, vertex_z);
437  return eta;
438 }
float getEta(const GlobalPoint &position, const float &vertex_z=0.) const
Definition: RecHitTools.cc:429
GlobalPoint getPosition(const DetId &id) const
Definition: RecHitTools.cc:129
static int position[264][3]
Definition: ReadPGInfo.cc:289
void RecHitTools::getEvent ( const edm::Event ev)

Definition at line 68 of file RecHitTools.cc.

Referenced by hgcal::ClusterTools::getEvent().

68 {}
void RecHitTools::getEventSetup ( const edm::EventSetup es)

Definition at line 70 of file RecHitTools.cc.

References bhLastLayer_, bhMaxIphi_, bhOffset_, fhLastLayer_, fhOffset_, DetId::Forward, ForwardEmpty, relativeConstraints::geom, geom_, geometryType_, edm::EventSetup::get(), CaloGeometry::getSubdetectorGeometry(), DetId::Hcal, HcalEndcap, HFNose, DetId::HGCalEE, DetId::HGCalHSc, DetId::HGCalHSi, HGCEE, HGCHEF, hgcalTopologyTester_cfi::layers, SiStripPI::max, maxNumberOfWafersNose_, maxNumberOfWafersPerLayer_, and edm::ESHandle< T >::product().

Referenced by HGCalHitCalibration::analyze(), HGCalShowerSeparation::analyze(), TICLLayerTileProducer::beginRun(), FilteredLayerClustersProducer::beginRun(), HGCalEgammaIDHelper::eventInit(), hgcal::ClusterTools::getEventSetup(), HGCalDepthPreClusterer::getEventSetup(), HGCal3DClustering::getEventSetup(), PFHGCalRecHitCreator< DET, Layer, det, subdet >::importRecHits(), ticl::PatternRecognitionbyCA::makeTracksters(), HGCalRecHitAbsAlgo::set(), and RealisticSimClusterMapper::update().

70  {
72  es.get<CaloGeometryRecord>().get(geom);
73 
74  geom_ = geom.product();
75  unsigned int wmaxEE(0), wmaxFH(0);
76  auto geomEE = static_cast<const HGCalGeometry*>(
78  //check if it's the new geometry
79  if (geomEE) {
80  geometryType_ = 1;
81  fhOffset_ = (geomEE->topology().dddConstants()).layers(true);
82  wmaxEE = (geomEE->topology().dddConstants()).waferCount(0);
83  auto geomFH = static_cast<const HGCalGeometry*>(
85  wmaxFH = (geomFH->topology().dddConstants()).waferCount(0);
86  fhLastLayer_ = fhOffset_ + (geomFH->topology().dddConstants()).layers(true);
87  auto geomBH = static_cast<const HGCalGeometry*>(
89  bhOffset_ =
90  fhOffset_ + (geomBH->topology().dddConstants()).firstLayer() - (geomEE->topology().dddConstants()).firstLayer();
91  bhLastLayer_ = bhOffset_ + (geomBH->topology().dddConstants()).layers(true);
92  bhMaxIphi_ = geomBH->topology().dddConstants().maxCells(true);
93  } else {
94  geometryType_ = 0;
95  geomEE =
97  fhOffset_ = (geomEE->topology().dddConstants()).layers(true);
98  wmaxEE = 1 + (geomEE->topology().dddConstants()).waferMax();
99  auto geomFH =
101  bhOffset_ = fhOffset_ + (geomFH->topology().dddConstants()).layers(true);
102  wmaxFH = 1 + (geomFH->topology().dddConstants()).waferMax();
104  auto geomBH =
106  bhLastLayer_ = bhOffset_ + (geomBH->topology().dddConstants())->getMaxDepth(1);
107  }
108  maxNumberOfWafersPerLayer_ = std::max(wmaxEE, wmaxFH);
109  // For nose geometry
110  auto geomNose =
112  if (geomNose) {
113  maxNumberOfWafersNose_ = (geomNose->topology().dddConstants()).waferCount(0);
114  } else {
116  }
117 }
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:34
unsigned int fhLastLayer_
Definition: RecHitTools.h:75
unsigned int bhLastLayer_
Definition: RecHitTools.h:75
unsigned int maxNumberOfWafersNose_
Definition: RecHitTools.h:76
unsigned int fhOffset_
Definition: RecHitTools.h:75
unsigned int bhOffset_
Definition: RecHitTools.h:75
unsigned int maxNumberOfWafersPerLayer_
Definition: RecHitTools.h:76
const CaloGeometry * geom_
Definition: RecHitTools.h:74
T get() const
Definition: EventSetup.h:73
T const * product() const
Definition: ESHandle.h:86
const CaloGeometry* hgcal::RecHitTools::getGeometry ( ) const
inline

Definition at line 61 of file RecHitTools.h.

61 { return geom_; };
const CaloGeometry * geom_
Definition: RecHitTools.h:74
int hgcal::RecHitTools::getGeometryType ( ) const
inline

Definition at line 70 of file RecHitTools.h.

Referenced by RealisticSimClusterMapper::buildClusters().

70 { return geometryType_; }
unsigned int RecHitTools::getLayer ( DetId::Detector  type,
bool  nose = false 
) const

Definition at line 297 of file RecHitTools.cc.

References DetId::Forward, ForwardEmpty, geom_, CaloGeometry::getSubdetectorGeometry(), HFNose, DetId::HGCalEE, DetId::HGCalHSc, DetId::HGCalHSi, and hgcalTopologyTester_cfi::layers.

Referenced by RealisticSimClusterMapper::buildClusters(), and getLayerWithOffset().

297  {
298  int layer;
299  switch (type) {
300  case (DetId::HGCalEE): {
301  auto geomEE =
303  layer = (geomEE->topology().dddConstants()).layers(true);
304  break;
305  }
306  case (DetId::HGCalHSi): {
307  auto geomFH =
309  layer = (geomFH->topology().dddConstants()).layers(true);
310  break;
311  }
312  case (DetId::HGCalHSc): {
313  auto geomBH =
315  layer = (geomBH->topology().dddConstants()).layers(true);
316  break;
317  }
318  case (DetId::Forward): {
319  if (nose) {
320  auto geomHFN = static_cast<const HGCalGeometry*>(
322  layer = (geomHFN->topology().dddConstants()).layers(true);
323  } else {
324  auto geomEE = static_cast<const HGCalGeometry*>(
326  layer = (geomEE->topology().dddConstants()).layers(true);
327  auto geomFH = static_cast<const HGCalGeometry*>(
329  layer += (geomFH->topology().dddConstants()).layers(true);
330  }
331  break;
332  }
333  default:
334  layer = 0;
335  }
336  return (unsigned int)(layer);
337 }
type
Definition: HCALResponse.h:21
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:34
const CaloGeometry * geom_
Definition: RecHitTools.h:74
unsigned int RecHitTools::getLayer ( ForwardSubdetector  type) const

Definition at line 247 of file RecHitTools.cc.

References DetId::Forward, ForwardEmpty, geom_, CaloGeometry::getSubdetectorGeometry(), DetId::Hcal, HcalEndcap, HFNose, HGCEE, HGCHEB, HGCHEF, and hgcalTopologyTester_cfi::layers.

247  {
248  int layer(0);
249  switch (type) {
250  case (ForwardSubdetector::HGCEE): {
251  auto geomEE =
253  layer = (geomEE->topology().dddConstants()).layers(true);
254  break;
255  }
257  auto geomFH =
259  layer = (geomFH->topology().dddConstants()).layers(true);
260  break;
261  }
263  auto geomBH =
265  layer = (geomBH->topology().dddConstants())->getMaxDepth(1);
266  break;
267  }
269  auto geomEE =
271  layer = (geomEE->topology().dddConstants()).layers(true);
272  auto geomFH =
274  layer += (geomFH->topology().dddConstants()).layers(true);
275  auto geomBH =
277  if (geomBH)
278  layer += (geomBH->topology().dddConstants())->getMaxDepth(1);
279  auto geomBH2 =
281  if (geomBH2)
282  layer += (geomBH2->topology().dddConstants()).layers(true);
283  break;
284  }
286  auto geomHFN =
288  layer = (geomHFN->topology().dddConstants()).layers(true);
289  break;
290  }
291  default:
292  layer = 0;
293  }
294  return (unsigned int)(layer);
295 }
type
Definition: HCALResponse.h:21
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:34
const CaloGeometry * geom_
Definition: RecHitTools.h:74
unsigned int RecHitTools::getLayer ( const DetId id) const

Definition at line 339 of file RecHitTools.cc.

References HcalDetId::depth(), DetId::Forward, DetId::Hcal, HcalEndcap, HFNose, DetId::HGCalEE, DetId::HGCalHSc, DetId::HGCalHSi, HGCalDetId::layer(), HGCScintillatorDetId::layer(), HFNoseDetId::layer(), HGCSiliconDetId::layer(), and SiStripPI::max.

339  {
340  unsigned int layer = std::numeric_limits<unsigned int>::max();
341  if (id.det() == DetId::HGCalEE || id.det() == DetId::HGCalHSi) {
342  layer = HGCSiliconDetId(id).layer();
343  } else if (id.det() == DetId::HGCalHSc) {
344  layer = HGCScintillatorDetId(id).layer();
345  } else if (id.det() == DetId::Forward && id.subdetId() == static_cast<int>(HFNose)) {
346  layer = HFNoseDetId(id).layer();
347  } else if (id.det() == DetId::Forward) {
348  layer = HGCalDetId(id).layer();
349  } else if (id.det() == DetId::Hcal && id.subdetId() == HcalEndcap) {
350  layer = HcalDetId(id).depth();
351  }
352  return layer;
353 }
int depth() const
get the tower depth
Definition: HcalDetId.h:164
int layer() const
get the layer #
Definition: HFNoseDetId.h:55
int layer() const
get the layer #
int layer() const
get the layer #
int layer() const
get the layer #
Definition: HGCalDetId.h:46
unsigned int RecHitTools::getLayerWithOffset ( const DetId id) const

Definition at line 355 of file RecHitTools.cc.

References bhOffset_, fhOffset_, DetId::Forward, getLayer(), DetId::Hcal, HcalEndcap, DetId::HGCalHSc, DetId::HGCalHSi, and HGCHEF.

Referenced by HGCalHitCalibration::analyze(), HGCalShowerSeparation::analyze(), RealisticSimClusterMapper::buildClusters(), ticl::PatternRecognitionbyCA::energyRegressionAndID(), HGCalHitCalibration::fillWithRecHits(), hgcal::ClusterTools::getLayer(), HGCalRecHitSimpleAlgo::makeRecHit(), HGCal3DClustering::organizeByLayer(), TICLLayerTileProducer::produce(), HGCalIsoCalculator::produceHGCalIso(), and hgcal::EGammaPCAHelper::storeRecHits().

355  {
356  unsigned int layer = getLayer(id);
357  if (id.det() == DetId::Forward && id.subdetId() == HGCHEF) {
358  layer += fhOffset_;
359  } else if (id.det() == DetId::HGCalHSi || id.det() == DetId::HGCalHSc) {
360  // DetId::HGCalHSc hits include the offset w.r.t. EE already
361  layer += fhOffset_;
362  } else if (id.det() == DetId::Hcal && id.subdetId() == HcalEndcap) {
363  layer += bhOffset_;
364  }
365  return layer;
366 }
unsigned int fhOffset_
Definition: RecHitTools.h:75
unsigned int bhOffset_
Definition: RecHitTools.h:75
unsigned int getLayer(DetId::Detector type, bool nose=false) const
Definition: RecHitTools.cc:297
float RecHitTools::getPhi ( const GlobalPoint position) const

Definition at line 440 of file RecHitTools.cc.

References PV3DBase< T, PVType, FrameType >::x(), and PV3DBase< T, PVType, FrameType >::y().

Referenced by HGCalIsoCalculator::setRecHits().

440  {
441  float phi = atan2(position.y(), position.x());
442  return phi;
443 }
T y() const
Definition: PV3DBase.h:60
T x() const
Definition: PV3DBase.h:59
float RecHitTools::getPhi ( const DetId id) const

Definition at line 445 of file RecHitTools.cc.

References getPosition(), position, PV3DBase< T, PVType, FrameType >::x(), and PV3DBase< T, PVType, FrameType >::y().

445  {
447  float phi = atan2(position.y(), position.x());
448  return phi;
449 }
T y() const
Definition: PV3DBase.h:60
GlobalPoint getPosition(const DetId &id) const
Definition: RecHitTools.cc:129
static int position[264][3]
Definition: ReadPGInfo.cc:289
T x() const
Definition: PV3DBase.h:59
GlobalPoint RecHitTools::getPosition ( const DetId id) const

Definition at line 129 of file RecHitTools.cc.

References relativeConstraints::geom, CaloGeometry::getGeometry(), HGCalGeometry::getPosition(), getSubdetectorGeometry(), DetId::Hcal, and position.

Referenced by HGCalShowerSeparation::analyze(), RealisticSimClusterMapper::buildClusters(), getEta(), getPhi(), getPt(), hgcal::ClusterTools::getWidths(), HGCalImagingAlgo::Hexel::Hexel(), HGCalIsoCalculator::setRecHits(), and hgcal::EGammaPCAHelper::storeRecHits().

129  {
130  auto geom = getSubdetectorGeometry(id);
132  if (id.det() == DetId::Hcal) {
133  position = geom->getGeometry(id)->getPosition();
134  } else {
135  auto hg = static_cast<const HGCalGeometry*>(geom);
136  position = hg->getPosition(id);
137  }
138  return position;
139 }
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
Definition: RecHitTools.cc:119
static int position[264][3]
Definition: ReadPGInfo.cc:289
GlobalPoint RecHitTools::getPositionLayer ( int  layer,
bool  nose = false 
) const

Definition at line 141 of file RecHitTools.cc.

References funct::abs(), fhOffset_, DetId::Forward, geom_, geometryType_, CaloGeometry::getSubdetectorGeometry(), HFNose, and HGCalDDDConstants::waferZ().

Referenced by hgcal::EGammaPCAHelper::findZFirstLayer().

141  {
142  unsigned int lay = std::abs(layer);
143  double z(0);
144  if (nose) {
145  auto geomNose =
147  if (geomNose) {
148  const HGCalDDDConstants* ddd = &(geomNose->topology().dddConstants());
149  if (ddd)
150  z = (layer > 0) ? ddd->waferZ(lay, true) : -ddd->waferZ(lay, true);
151  }
152  } else {
153  const HGCalDDDConstants* ddd = get_ddd(geom_, geometryType_, fhOffset_, lay);
154  if (geometryType_ == 1) {
155  if (lay > fhOffset_)
156  lay -= fhOffset_;
157  }
158  z = (layer > 0) ? ddd->waferZ(lay, true) : -ddd->waferZ(lay, true);
159  }
160  return GlobalPoint(0, 0, z);
161 }
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:34
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
unsigned int fhOffset_
Definition: RecHitTools.h:75
const CaloGeometry * geom_
Definition: RecHitTools.h:74
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double waferZ(int layer, bool reco) const
float RecHitTools::getPt ( const GlobalPoint position,
const float &  hitEnergy,
const float &  vertex_z = 0. 
) const

Definition at line 451 of file RecHitTools.cc.

References PVValHelper::eta, getEta(), and DiDispStaMuonMonitor_cfi::pt.

451  {
452  float eta = getEta(position, vertex_z);
453  float pt = hitEnergy / cosh(eta);
454  return pt;
455 }
float getEta(const GlobalPoint &position, const float &vertex_z=0.) const
Definition: RecHitTools.cc:429
float RecHitTools::getPt ( const DetId id,
const float &  hitEnergy,
const float &  vertex_z = 0. 
) const

Definition at line 457 of file RecHitTools.cc.

References PVValHelper::eta, getEta(), getPosition(), position, and DiDispStaMuonMonitor_cfi::pt.

457  {
459  float eta = getEta(position, vertex_z);
460  float pt = hitEnergy / cosh(eta);
461  return pt;
462 }
float getEta(const GlobalPoint &position, const float &vertex_z=0.) const
Definition: RecHitTools.cc:429
GlobalPoint getPosition(const DetId &id) const
Definition: RecHitTools.cc:129
static int position[264][3]
Definition: ReadPGInfo.cc:289
std::float_t RecHitTools::getRadiusToSide ( const DetId id) const

Definition at line 225 of file RecHitTools.cc.

References HGCalDDDConstants::cellSizeHex(), TauDecayModes::dec, DetId::Forward, getSubdetectorGeometry(), HFNose, DetId::HGCalEE, DetId::HGCalHSi, SiStripPI::max, findQualityFiles::size, HFNoseDetId::type(), HGCSiliconDetId::type(), and HGCalDetId::waferType().

225  {
226  auto geom = getSubdetectorGeometry(id);
228  if (id.det() == DetId::HGCalEE || id.det() == DetId::HGCalHSi) {
229  const HGCSiliconDetId hid(id);
230  auto ddd = get_ddd(geom, hid);
231  size = ddd->cellSizeHex(hid.type());
232  } else if (id.det() == DetId::Forward && id.subdetId() == static_cast<int>(HFNose)) {
233  const HFNoseDetId hid(id);
234  auto ddd = get_ddd(geom, hid);
235  size = ddd->cellSizeHex(hid.type());
236  } else if (id.det() == DetId::Forward) {
237  const HGCalDetId hid(id);
238  auto ddd = get_ddd(geom, hid);
239  size = ddd->cellSizeHex(hid.waferType());
240  } else {
241  edm::LogError("getRadiusToSide::InvalidSiliconDetid")
242  << "det id: " << std::hex << id.rawId() << std::dec << ":" << id.det() << " is not HGCal silicon!";
243  }
244  return size;
245 }
size
Write out results.
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
Definition: RecHitTools.cc:119
int hgcal::RecHitTools::getScintMaxIphi ( ) const
inline

Definition at line 69 of file RecHitTools.h.

69 { return bhMaxIphi_; }
int RecHitTools::getSiThickIndex ( const DetId id) const

Definition at line 205 of file RecHitTools.cc.

References DetId::Forward, getSiThickness(), HFNose, DetId::HGCalEE, DetId::HGCalHSi, Calorimetry_cff::thickness, HFNoseDetId::type(), and HGCSiliconDetId::type().

Referenced by hgcal::EGammaPCAHelper::storeRecHits().

205  {
206  int thickIndex = -1;
207  if (id.det() == DetId::HGCalEE || id.det() == DetId::HGCalHSi) {
208  thickIndex = HGCSiliconDetId(id).type();
209  } else if (id.det() == DetId::Forward && id.subdetId() == static_cast<int>(HFNose)) {
210  thickIndex = HFNoseDetId(id).type();
211  } else if (id.det() == DetId::Forward) {
212  float thickness = getSiThickness(id);
213  if (thickness > 99. && thickness < 121.)
214  thickIndex = 0;
215  else if (thickness > 199. && thickness < 201.)
216  thickIndex = 1;
217  else if (thickness > 299. && thickness < 301.)
218  thickIndex = 2;
219  else
220  assert(thickIndex > 0 && "ERROR - silicon thickness has a nonsensical value");
221  }
222  return thickIndex;
223 }
int type() const
get the type
Definition: HFNoseDetId.h:49
int type() const
get the type
std::float_t getSiThickness(const DetId &) const
Definition: RecHitTools.cc:179
std::float_t RecHitTools::getSiThickness ( const DetId id) const

Definition at line 179 of file RecHitTools.cc.

References HGCalDDDConstants::cellThickness(), TauDecayModes::dec, DetId::Forward, getSubdetectorGeometry(), HFNose, DetId::HGCalEE, DetId::HGCalHSi, HGCalDetId::layer(), HFNoseDetId::layer(), HGCSiliconDetId::layer(), LogDebug, HGCalDetId::wafer(), HFNoseDetId::waferU(), HGCSiliconDetId::waferU(), HFNoseDetId::waferV(), and HGCSiliconDetId::waferV().

Referenced by HGCalHitCalibration::analyze(), HGCalHitCalibration::fillWithRecHits(), and getSiThickIndex().

179  {
180  auto geom = getSubdetectorGeometry(id);
181  std::float_t thick(0.37);
182  if (id.det() == DetId::HGCalEE || id.det() == DetId::HGCalHSi) {
183  const HGCSiliconDetId hid(id);
184  auto ddd = get_ddd(geom, hid);
185  thick = ddd->cellThickness(hid.layer(), hid.waferU(), hid.waferV());
186  } else if (id.det() == DetId::Forward && id.subdetId() == static_cast<int>(HFNose)) {
187  const HFNoseDetId hid(id);
188  auto ddd = get_ddd(geom, hid);
189  thick = ddd->cellThickness(hid.layer(), hid.waferU(), hid.waferV());
190  } else if (id.det() == DetId::Forward) {
191  const HGCalDetId hid(id);
192  auto ddd = get_ddd(geom, hid);
193  thick = ddd->cellThickness(hid.layer(), hid.wafer(), 0);
194  } else {
195  LogDebug("getSiThickness::InvalidSiliconDetid")
196  << "det id: " << std::hex << id.rawId() << std::dec << ":" << id.det() << " is not HGCal silicon!";
197  // It does not make any sense to return 0.37 as thickness for a detector
198  // that is not Silicon(does it?). Mark it as 0. to be able to distinguish
199  // the case.
200  thick = 0.f;
201  }
202  return thick;
203 }
#define LogDebug(id)
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
Definition: RecHitTools.cc:119
const CaloSubdetectorGeometry * RecHitTools::getSubdetectorGeometry ( const DetId id) const

Definition at line 119 of file RecHitTools.cc.

References ForwardEmpty, relativeConstraints::geom, geom_, CaloGeometry::getSubdetectorGeometry(), DetId::HGCalEE, DetId::HGCalHSc, and DetId::HGCalHSi.

Referenced by getPosition(), getRadiusToSide(), getSiThickness(), isHalfCell(), and maskCell().

119  {
120  DetId::Detector det = id.det();
121  int subdet = (det == DetId::HGCalEE || det == DetId::HGCalHSi || det == DetId::HGCalHSc)
123  : id.subdetId();
124  auto geom = geom_->getSubdetectorGeometry(det, subdet);
125  check_geom(geom);
126  return geom;
127 }
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:34
const CaloGeometry * geom_
Definition: RecHitTools.h:74
Detector
Definition: DetId.h:24
std::pair< int, int > RecHitTools::getWafer ( const DetId id) const

Definition at line 368 of file RecHitTools.cc.

References TauDecayModes::dec, DetId::Forward, HFNose, DetId::HGCalEE, DetId::HGCalHSi, SiStripPI::max, HGCalDetId::wafer(), HFNoseDetId::waferU(), HGCSiliconDetId::waferU(), HFNoseDetId::waferV(), and HGCSiliconDetId::waferV().

368  {
369  int waferU = std::numeric_limits<int>::max();
370  int waferV = 0;
371  if ((id.det() == DetId::HGCalEE) || (id.det() == DetId::HGCalHSi)) {
372  waferU = HGCSiliconDetId(id).waferU();
373  waferV = HGCSiliconDetId(id).waferV();
374  } else if (id.det() == DetId::Forward && id.subdetId() == static_cast<int>(HFNose)) {
375  waferU = HFNoseDetId(id).waferU();
376  waferV = HFNoseDetId(id).waferV();
377  } else if (id.det() == DetId::Forward) {
378  waferU = HGCalDetId(id).wafer();
379  } else {
380  edm::LogError("getWafer::InvalidSiliconDetid")
381  << "det id: " << std::hex << id.rawId() << std::dec << ":" << id.det() << " is not HGCal silicon!";
382  }
383  return std::pair<int, int>(waferU, waferV);
384 }
int waferU() const
int waferV() const
Definition: HFNoseDetId.h:77
int waferU() const
Definition: HFNoseDetId.h:74
int waferV() const
int wafer() const
get the wafer #
Definition: HGCalDetId.h:40
bool RecHitTools::isHalfCell ( const DetId id) const

Definition at line 404 of file RecHitTools.cc.

References HGCalDetId::cell(), DetId::Forward, getSubdetectorGeometry(), HGCalDDDConstants::isHalfCell(), HGCalDetId::waferType(), and HGCalDDDConstants::waferTypeT().

404  {
405  bool ishalf = false;
406  if (id.det() == DetId::Forward) {
407  HGCalDetId hid(id);
408  auto geom = getSubdetectorGeometry(hid);
409  auto ddd = get_ddd(geom, hid);
410  const int waferType = ddd->waferTypeT(hid.waferType());
411  return ddd->isHalfCell(waferType, hid.cell());
412  }
413  //new geometry is always false
414  return ishalf;
415 }
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
Definition: RecHitTools.cc:119
bool RecHitTools::isOnlySilicon ( const unsigned int  layer) const

Definition at line 424 of file RecHitTools.cc.

References bhLastLayer_, and bhOffset_.

424  {
425  bool isonlysilicon = (layer % bhLastLayer_) < bhOffset_;
426  return isonlysilicon;
427 }
unsigned int bhLastLayer_
Definition: RecHitTools.h:75
unsigned int bhOffset_
Definition: RecHitTools.h:75
bool RecHitTools::isSilicon ( const DetId id) const

Definition at line 417 of file RecHitTools.cc.

References DetId::HGCalEE, and DetId::HGCalHSi.

Referenced by ticl::ClusterFilterByAlgoAndSize::filter().

417  {
418  bool issilicon = false;
419  if (id.det() == DetId::HGCalEE || id.det() == DetId::HGCalHSi)
420  issilicon = true;
421  return issilicon;
422 }
unsigned int hgcal::RecHitTools::lastLayerBH ( ) const
inline

Definition at line 65 of file RecHitTools.h.

Referenced by HGCal3DClustering::getEventSetup(), and hgcal::EGammaPCAHelper::setRecHitTools().

65 { return bhLastLayer_; }
unsigned int bhLastLayer_
Definition: RecHitTools.h:75
unsigned int hgcal::RecHitTools::lastLayerEE ( ) const
inline
unsigned int hgcal::RecHitTools::lastLayerFH ( ) const
inline
bool RecHitTools::maskCell ( const DetId id,
int  corners = 3 
) const

Definition at line 464 of file RecHitTools.cc.

References HGCalTopology::dddConstants(), getSubdetectorGeometry(), DetId::Hcal, HGCalDDDConstants::maskCell(), and HGCalGeometry::topology().

464  {
465  if (id.det() == DetId::Hcal) {
466  return false;
467  } else {
468  auto hg = static_cast<const HGCalGeometry*>(getSubdetectorGeometry(id));
469  return hg->topology().dddConstants().maskCell(id, corners);
470  }
471 }
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
Definition: RecHitTools.cc:119
const HGCalTopology & topology() const
const HGCalDDDConstants & dddConstants() const
Definition: HGCalTopology.h:96
bool maskCell(const DetId &id, int corners) const
unsigned int hgcal::RecHitTools::maxNumberOfWafersPerLayer ( bool  nose = false) const
inline

Definition at line 66 of file RecHitTools.h.

66  {
68  }
unsigned int maxNumberOfWafersNose_
Definition: RecHitTools.h:76
unsigned int maxNumberOfWafersPerLayer_
Definition: RecHitTools.h:76
int RecHitTools::zside ( const DetId id) const

Definition at line 163 of file RecHitTools.cc.

References DetId::Forward, DetId::Hcal, HcalEndcap, HFNose, DetId::HGCalEE, DetId::HGCalHSc, DetId::HGCalHSi, HGCScintillatorDetId::zside(), HGCalDetId::zside(), HFNoseDetId::zside(), HGCSiliconDetId::zside(), and HcalDetId::zside().

Referenced by TICLLayerTileProducer::produce().

163  {
164  int zside = 0;
165  if (id.det() == DetId::HGCalEE || id.det() == DetId::HGCalHSi) {
166  zside = HGCSiliconDetId(id).zside();
167  } else if (id.det() == DetId::HGCalHSc) {
168  zside = HGCScintillatorDetId(id).zside();
169  } else if (id.det() == DetId::Forward && id.subdetId() == static_cast<int>(HFNose)) {
170  zside = HFNoseDetId(id).zside();
171  } else if (id.det() == DetId::Forward) {
172  zside = HGCalDetId(id).zside();
173  } else if (id.det() == DetId::Hcal && id.subdetId() == HcalEndcap) {
174  zside = HcalDetId(id).zside();
175  }
176  return zside;
177 }
int zside() const
get the z-side of the cell (1/-1)
Definition: HFNoseDetId.h:52
int zside(const DetId &id) const
Definition: RecHitTools.cc:163
int zside() const
get the z-side of the cell (1/-1)
Definition: HcalDetId.h:141
int zside() const
get the z-side of the cell (1/-1)
int zside() const
get the z-side of the cell (1/-1)
Definition: HGCalDetId.h:49
int zside() const
get the z-side of the cell (1/-1)

Member Data Documentation

unsigned int hgcal::RecHitTools::bhLastLayer_
private

Definition at line 75 of file RecHitTools.h.

Referenced by getEventSetup(), and isOnlySilicon().

int hgcal::RecHitTools::bhMaxIphi_
private

Definition at line 78 of file RecHitTools.h.

Referenced by getEventSetup().

unsigned int hgcal::RecHitTools::bhOffset_
private

Definition at line 75 of file RecHitTools.h.

Referenced by getEventSetup(), getLayerWithOffset(), and isOnlySilicon().

unsigned int hgcal::RecHitTools::fhLastLayer_
private

Definition at line 75 of file RecHitTools.h.

Referenced by getEventSetup().

unsigned int hgcal::RecHitTools::fhOffset_
private

Definition at line 75 of file RecHitTools.h.

Referenced by getEventSetup(), getLayerWithOffset(), and getPositionLayer().

const CaloGeometry* hgcal::RecHitTools::geom_
private

Definition at line 74 of file RecHitTools.h.

Referenced by getEventSetup(), getLayer(), getPositionLayer(), and getSubdetectorGeometry().

int hgcal::RecHitTools::geometryType_
private

Definition at line 77 of file RecHitTools.h.

Referenced by getEventSetup(), and getPositionLayer().

unsigned int hgcal::RecHitTools::maxNumberOfWafersNose_
private

Definition at line 76 of file RecHitTools.h.

Referenced by getEventSetup().

unsigned int hgcal::RecHitTools::maxNumberOfWafersPerLayer_
private

Definition at line 76 of file RecHitTools.h.

Referenced by getEventSetup().