CMS 3D CMS Logo

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

#include <RecHitTools.h>

Public Member Functions

std::pair< uint32_t, uint32_t > firstAndLastLayer (DetId::Detector det, int subdet) const
 
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
 
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
 
std::pair< float, float > getScintDEtaDPhi (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 isScintillator (const DetId &) const
 
bool isSilicon (const DetId &) const
 
unsigned int lastLayer (bool nose=false) const
 
unsigned int lastLayerBH () const
 
unsigned int lastLayerEE (bool nose=false) const
 
unsigned int lastLayerFH () const
 
bool maskCell (const DetId &id, int corners=3) const
 
unsigned int maxNumberOfWafersPerLayer (bool nose=false) const
 
 RecHitTools ()
 
void setGeometry (CaloGeometry const &)
 
int zside (const DetId &id) const
 
 ~RecHitTools ()
 

Private Attributes

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

Detailed Description

Definition at line 23 of file RecHitTools.h.

Constructor & Destructor Documentation

◆ RecHitTools()

hgcal::RecHitTools::RecHitTools ( )
inline

Definition at line 25 of file RecHitTools.h.

26  : geom_(nullptr),
27  eeOffset_(0),
28  fhOffset_(0),
29  bhFirstLayer_(0),
30  bhOffset_(0),
31  fhLastLayer_(0),
32  noseLastLayer_(0),
33  geometryType_(0) {}
unsigned int bhFirstLayer_
Definition: RecHitTools.h:90
unsigned int fhLastLayer_
Definition: RecHitTools.h:90
unsigned int fhOffset_
Definition: RecHitTools.h:90
unsigned int bhOffset_
Definition: RecHitTools.h:90
const CaloGeometry * geom_
Definition: RecHitTools.h:89
unsigned int eeOffset_
Definition: RecHitTools.h:90
unsigned int noseLastLayer_
Definition: RecHitTools.h:90

◆ ~RecHitTools()

hgcal::RecHitTools::~RecHitTools ( )
inline

Definition at line 34 of file RecHitTools.h.

34 {}

Member Function Documentation

◆ firstAndLastLayer()

std::pair< uint32_t, uint32_t > RecHitTools::firstAndLastLayer ( DetId::Detector  det,
int  subdet 
) const

Definition at line 476 of file RecHitTools.cc.

References bhFirstLayer_, bhLastLayer_, eeOffset_, fhLastLayer_, fhOffset_, DetId::Forward, HFNose, DetId::HGCalEE, DetId::HGCalHSi, HGCEE, HGCHEF, and noseLastLayer_.

476  {
477  if ((det == DetId::HGCalEE) || ((det == DetId::Forward) && (subdet == HGCEE))) {
478  return std::make_pair(eeOffset_ + 1, fhOffset_);
479  } else if ((det == DetId::HGCalHSi) || ((det == DetId::Forward) && (subdet == HGCHEF))) {
480  return std::make_pair(fhOffset_ + 1, fhLastLayer_);
481  } else if ((det == DetId::Forward) && (subdet == HFNose)) {
482  return std::make_pair(1, noseLastLayer_);
483  } else {
484  return std::make_pair(bhFirstLayer_, bhLastLayer_);
485  }
486 }
unsigned int bhFirstLayer_
Definition: RecHitTools.h:90
unsigned int fhLastLayer_
Definition: RecHitTools.h:90
unsigned int bhLastLayer_
Definition: RecHitTools.h:90
unsigned int fhOffset_
Definition: RecHitTools.h:90
unsigned int eeOffset_
Definition: RecHitTools.h:90
unsigned int noseLastLayer_
Definition: RecHitTools.h:90

◆ firstLayerBH()

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

◆ getCell()

std::pair< int, int > RecHitTools::getCell ( const DetId id) const

Definition at line 397 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.

397  {
398  int cellU = std::numeric_limits<int>::max();
399  int cellV = 0;
400  if ((id.det() == DetId::HGCalEE) || (id.det() == DetId::HGCalHSi)) {
401  cellU = HGCSiliconDetId(id).cellU();
402  cellV = HGCSiliconDetId(id).cellV();
403  } else if (id.det() == DetId::Forward && id.subdetId() == static_cast<int>(HFNose)) {
404  cellU = HFNoseDetId(id).cellU();
405  cellV = HFNoseDetId(id).cellV();
406  } else if (id.det() == DetId::Forward) {
407  cellU = HGCalDetId(id).cell();
408  } else {
409  edm::LogError("getCell::InvalidSiliconDetid")
410  << "det id: " << std::hex << id.rawId() << std::dec << ":" << id.det() << " is not HGCal silicon!";
411  }
412  return std::pair<int, int>(cellU, cellV);
413 }
int cellU() const
get the cell #&#39;s in u,v or in x,y
Log< level::Error, false > LogError
int cellU() const
get the cell #&#39;s in u,v or in x,y
Definition: HFNoseDetId.h:59
int cell() const
get the absolute value of the cell #&#39;s in x and y
Definition: HGCalDetId.h:37
int cellV() const
Definition: HFNoseDetId.h:60
int cellV() const

◆ getEta() [1/2]

float RecHitTools::getEta ( const GlobalPoint position,
const float &  vertex_z = 0. 
) const

Definition at line 441 of file RecHitTools.cc.

References PV3DBase< T, PVType, FrameType >::eta(), and position.

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

441  {
442  GlobalPoint corrected_position = GlobalPoint(position.x(), position.y(), position.z() - vertex_z);
443  return corrected_position.eta();
444 }
T eta() const
Definition: PV3DBase.h:73
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
static int position[264][3]
Definition: ReadPGInfo.cc:289

◆ getEta() [2/2]

float RecHitTools::getEta ( const DetId id,
const float &  vertex_z = 0. 
) const

Definition at line 446 of file RecHitTools.cc.

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

446  {
448  float eta = getEta(position, vertex_z);
449  return eta;
450 }
GlobalPoint getPosition(const DetId &id) const
Definition: RecHitTools.cc:129
float getEta(const GlobalPoint &position, const float &vertex_z=0.) const
Definition: RecHitTools.cc:441
static int position[264][3]
Definition: ReadPGInfo.cc:289

◆ getGeometry()

const CaloGeometry* hgcal::RecHitTools::getGeometry ( ) const
inline

Definition at line 74 of file RecHitTools.h.

References geom_.

74 { return geom_; };
const CaloGeometry * geom_
Definition: RecHitTools.h:89

◆ getGeometryType()

int hgcal::RecHitTools::getGeometryType ( ) const
inline

Definition at line 85 of file RecHitTools.h.

References geometryType_.

Referenced by RealisticSimClusterMapper::buildClusters().

85 { return geometryType_; }

◆ getLayer() [1/3]

unsigned int RecHitTools::getLayer ( DetId::Detector  type,
bool  nose = false 
) const

Definition at line 307 of file RecHitTools.cc.

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

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

307  {
308  int layer;
309  switch (type) {
310  case (DetId::HGCalEE): {
311  auto geomEE =
313  layer = (geomEE->topology().dddConstants()).layers(true);
314  break;
315  }
316  case (DetId::HGCalHSi): {
317  auto geomFH =
319  layer = (geomFH->topology().dddConstants()).layers(true);
320  break;
321  }
322  case (DetId::HGCalHSc): {
323  auto geomBH =
325  layer = (geomBH->topology().dddConstants()).layers(true);
326  break;
327  }
328  case (DetId::Forward): {
329  if (nose) {
330  auto geomHFN = static_cast<const HGCalGeometry*>(
332  layer = (geomHFN->topology().dddConstants()).layers(true);
333  } else {
334  auto geomEE = static_cast<const HGCalGeometry*>(
336  layer = (geomEE->topology().dddConstants()).layers(true);
337  auto geomFH = static_cast<const HGCalGeometry*>(
339  layer += (geomFH->topology().dddConstants()).layers(true);
340  }
341  break;
342  }
343  default:
344  layer = 0;
345  }
346  return (unsigned int)(layer);
347 }
constexpr std::array< uint8_t, layerIndexSize > layer
const CaloGeometry * geom_
Definition: RecHitTools.h:89
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:34

◆ getLayer() [2/3]

unsigned int RecHitTools::getLayer ( ForwardSubdetector  type) const

Definition at line 257 of file RecHitTools.cc.

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

257  {
258  int layer(0);
259  switch (type) {
260  case (ForwardSubdetector::HGCEE): {
261  auto geomEE =
263  layer = (geomEE->topology().dddConstants()).layers(true);
264  break;
265  }
267  auto geomFH =
269  layer = (geomFH->topology().dddConstants()).layers(true);
270  break;
271  }
273  auto geomBH =
275  layer = (geomBH->topology().dddConstants())->getMaxDepth(1);
276  break;
277  }
279  auto geomEE =
281  layer = (geomEE->topology().dddConstants()).layers(true);
282  auto geomFH =
284  layer += (geomFH->topology().dddConstants()).layers(true);
285  auto geomBH =
287  if (geomBH)
288  layer += (geomBH->topology().dddConstants())->getMaxDepth(1);
289  auto geomBH2 =
291  if (geomBH2)
292  layer += (geomBH2->topology().dddConstants()).layers(true);
293  break;
294  }
296  auto geomHFN =
298  layer = (geomHFN->topology().dddConstants()).layers(true);
299  break;
300  }
301  default:
302  layer = 0;
303  }
304  return (unsigned int)(layer);
305 }
constexpr std::array< uint8_t, layerIndexSize > layer
const CaloGeometry * geom_
Definition: RecHitTools.h:89
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:34

◆ getLayer() [3/3]

unsigned int RecHitTools::getLayer ( const DetId id) const

Definition at line 349 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(), phase1PixelTopology::layer, and SiStripPI::max.

349  {
351  if (id.det() == DetId::HGCalEE || id.det() == DetId::HGCalHSi) {
352  layer = HGCSiliconDetId(id).layer();
353  } else if (id.det() == DetId::HGCalHSc) {
355  } else if (id.det() == DetId::Forward && id.subdetId() == static_cast<int>(HFNose)) {
356  layer = HFNoseDetId(id).layer();
357  } else if (id.det() == DetId::Forward) {
358  layer = HGCalDetId(id).layer();
359  } else if (id.det() == DetId::Hcal && id.subdetId() == HcalEndcap) {
360  layer = HcalDetId(id).depth();
361  }
362  return layer;
363 }
int layer() const
get the layer #
Definition: HFNoseDetId.h:56
constexpr std::array< uint8_t, layerIndexSize > layer
int layer() const
get the layer #
int layer() const
get the layer #
int layer() const
get the layer #
Definition: HGCalDetId.h:46
constexpr int depth() const
get the tower depth
Definition: HcalDetId.h:164

◆ getLayerWithOffset()

unsigned int RecHitTools::getLayerWithOffset ( const DetId id) const

Definition at line 365 of file RecHitTools.cc.

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

Referenced by FWRecoGeometryESProducer::addCaloGeometry(), HGCalShowerSeparation::analyze(), HGCalHitCalibration::analyze(), RealisticSimClusterMapper::buildClusters(), TICLTrackstersEdgesValidation::dqmAnalyze(), TrackstersMergeProducer::energyRegressionAndID(), HGCalHitCalibration::fillWithRecHits(), ticl::ClusterFilterByAlgoAndSizeAndLayerRange::filter(), hgcal::ClusterTools::getLayer(), HGCalRecHitSimpleAlgo::makeRecHit(), HGCal3DClustering::organizeByLayer(), TICLLayerTileProducer::produce(), HGCalIsoCalculator::produceHGCalIso(), and hgcal::EGammaPCAHelper::storeRecHits().

365  {
366  unsigned int layer = getLayer(id);
367  if (id.det() == DetId::Forward && id.subdetId() == HGCHEF) {
368  layer += fhOffset_;
369  } else if (id.det() == DetId::HGCalHSi || id.det() == DetId::HGCalHSc) {
370  // DetId::HGCalHSc hits include the offset w.r.t. EE already
371  layer += fhOffset_;
372  } else if (id.det() == DetId::Hcal && id.subdetId() == HcalEndcap) {
373  layer += bhOffset_;
374  }
375  // no need to add offset for HFnose
376  return layer;
377 }
unsigned int getLayer(DetId::Detector type, bool nose=false) const
Definition: RecHitTools.cc:307
unsigned int fhOffset_
Definition: RecHitTools.h:90
unsigned int bhOffset_
Definition: RecHitTools.h:90
constexpr std::array< uint8_t, layerIndexSize > layer

◆ getPhi() [1/2]

float RecHitTools::getPhi ( const GlobalPoint position) const

Definition at line 452 of file RecHitTools.cc.

References position.

Referenced by HGCalIsoCalculator::setRecHits().

452  {
453  float phi = atan2(position.y(), position.x());
454  return phi;
455 }
static int position[264][3]
Definition: ReadPGInfo.cc:289

◆ getPhi() [2/2]

float RecHitTools::getPhi ( const DetId id) const

Definition at line 457 of file RecHitTools.cc.

References getPosition(), and position.

457  {
459  float phi = atan2(position.y(), position.x());
460  return phi;
461 }
GlobalPoint getPosition(const DetId &id) const
Definition: RecHitTools.cc:129
static int position[264][3]
Definition: ReadPGInfo.cc:289

◆ getPosition()

GlobalPoint RecHitTools::getPosition ( const DetId id) const

Definition at line 129 of file RecHitTools.cc.

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

Referenced by HGCalShowerSeparation::analyze(), RealisticSimClusterMapper::buildClusters(), getEta(), getPhi(), getPt(), hgcal::ClusterTools::getWidths(), 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

◆ getPositionLayer()

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, phase1PixelTopology::layer, and HGCalDDDConstants::waferZ().

Referenced by PositionAtECalEntranceComputer::beginEvent(), hgcal::EGammaPCAHelper::findZFirstLayer(), TrackstersMergeProducer::produce(), and SimTrackstersProducer::produce().

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 }
double waferZ(int layer, bool reco) const
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
unsigned int fhOffset_
Definition: RecHitTools.h:90
constexpr std::array< uint8_t, layerIndexSize > layer
const CaloGeometry * geom_
Definition: RecHitTools.h:89
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:34

◆ getPt() [1/2]

float RecHitTools::getPt ( const GlobalPoint position,
const float &  hitEnergy,
const float &  vertex_z = 0. 
) const

Definition at line 463 of file RecHitTools.cc.

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

463  {
464  float eta = getEta(position, vertex_z);
465  float pt = hitEnergy / cosh(eta);
466  return pt;
467 }
float getEta(const GlobalPoint &position, const float &vertex_z=0.) const
Definition: RecHitTools.cc:441
static int position[264][3]
Definition: ReadPGInfo.cc:289

◆ getPt() [2/2]

float RecHitTools::getPt ( const DetId id,
const float &  hitEnergy,
const float &  vertex_z = 0. 
) const

Definition at line 469 of file RecHitTools.cc.

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

469  {
471  float eta = getEta(position, vertex_z);
472  float pt = hitEnergy / cosh(eta);
473  return pt;
474 }
GlobalPoint getPosition(const DetId &id) const
Definition: RecHitTools.cc:129
float getEta(const GlobalPoint &position, const float &vertex_z=0.) const
Definition: RecHitTools.cc:441
static int position[264][3]
Definition: ReadPGInfo.cc:289

◆ getRadiusToSide()

std::float_t RecHitTools::getRadiusToSide ( const DetId id) const

Definition at line 235 of file RecHitTools.cc.

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

235  {
236  auto geom = getSubdetectorGeometry(id);
238  if (id.det() == DetId::HGCalEE || id.det() == DetId::HGCalHSi) {
239  const HGCSiliconDetId hid(id);
240  auto ddd = get_ddd(geom, hid);
241  size = ddd->cellSizeHex(hid.type());
242  } else if (id.det() == DetId::Forward && id.subdetId() == static_cast<int>(HFNose)) {
243  const HFNoseDetId hid(id);
244  auto ddd = get_ddd(geom, hid);
245  size = ddd->cellSizeHex(hid.type());
246  } else if (id.det() == DetId::Forward) {
247  const HGCalDetId hid(id);
248  auto ddd = get_ddd(geom, hid);
249  size = ddd->cellSizeHex(hid.waferType());
250  } else {
251  edm::LogError("getRadiusToSide::InvalidSiliconDetid")
252  << "det id: " << std::hex << id.rawId() << std::dec << ":" << id.det() << " is not HGCal silicon!";
253  }
254  return size;
255 }
size
Write out results.
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
Definition: RecHitTools.cc:119
Log< level::Error, false > LogError

◆ getScintDEtaDPhi()

std::pair< float, float > RecHitTools::getScintDEtaDPhi ( const DetId id) const

Definition at line 225 of file RecHitTools.cc.

References TauDecayModes::dec, CaloSubdetectorGeometry::getGeometry(), getSubdetectorGeometry(), isScintillator(), and LogDebug.

225  {
226  if (!isScintillator(id)) {
227  LogDebug("getScintDEtaDPhi::InvalidScintDetid")
228  << "det id: " << std::hex << id.rawId() << std::dec << ":" << id.det() << " is not HGCal scintillator!";
229  return {0.f, 0.f};
230  }
231  auto cellGeom = getSubdetectorGeometry(id)->getGeometry(id);
232  return {cellGeom->etaSpan(), cellGeom->phiSpan()};
233 }
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
Definition: RecHitTools.cc:119
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.
bool isScintillator(const DetId &) const
Definition: RecHitTools.cc:433
#define LogDebug(id)

◆ getScintMaxIphi()

int hgcal::RecHitTools::getScintMaxIphi ( ) const
inline

Definition at line 84 of file RecHitTools.h.

References bhMaxIphi_.

Referenced by HGCalClusteringAlgoBase::getEventSetup().

84 { return bhMaxIphi_; }

◆ getSiThickIndex()

int RecHitTools::getSiThickIndex ( const DetId id) const

Definition at line 205 of file RecHitTools.cc.

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

Referenced by FWRecoGeometryESProducer::addCaloGeometry(), and 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 }
assert(be >=bs)
int type() const
get the type
Definition: HFNoseDetId.h:50
int type() const
get the type
std::float_t getSiThickness(const DetId &) const
Definition: RecHitTools.cc:179

◆ getSiThickness()

std::float_t RecHitTools::getSiThickness ( const DetId id) const

Definition at line 179 of file RecHitTools.cc.

References HGCalDDDConstants::cellThickness(), TauDecayModes::dec, DetId::Forward, relativeConstraints::geom, getSubdetectorGeometry(), HFNose, DetId::HGCalEE, DetId::HGCalHSi, HGCalDetId::layer(), HFNoseDetId::layer(), HGCSiliconDetId::layer(), LogDebug, HGCalDetId::wafer(), HFNoseDetId::waferU(), HGCSiliconDetId::waferU(), HGCSiliconDetId::waferV(), and HFNoseDetId::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 }
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
Definition: RecHitTools.cc:119
#define LogDebug(id)

◆ getSubdetectorGeometry()

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(), getScintDEtaDPhi(), 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 CaloGeometry * geom_
Definition: RecHitTools.h:89
Detector
Definition: DetId.h:24
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:34

◆ getWafer()

std::pair< int, int > RecHitTools::getWafer ( const DetId id) const

Definition at line 379 of file RecHitTools.cc.

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

379  {
381  int waferV = 0;
382  if ((id.det() == DetId::HGCalEE) || (id.det() == DetId::HGCalHSi)) {
383  waferU = HGCSiliconDetId(id).waferU();
384  waferV = HGCSiliconDetId(id).waferV();
385  } else if (id.det() == DetId::Forward && id.subdetId() == static_cast<int>(HFNose)) {
386  waferU = HFNoseDetId(id).waferU();
387  waferV = HFNoseDetId(id).waferV();
388  } else if (id.det() == DetId::Forward) {
389  waferU = HGCalDetId(id).wafer();
390  } else {
391  edm::LogError("getWafer::InvalidSiliconDetid")
392  << "det id: " << std::hex << id.rawId() << std::dec << ":" << id.det() << " is not HGCal silicon!";
393  }
394  return std::pair<int, int>(waferU, waferV);
395 }
int wafer() const
get the wafer #
Definition: HGCalDetId.h:40
int32_t waferU(const int32_t index)
int waferU() const
Definition: HFNoseDetId.h:75
Log< level::Error, false > LogError
int waferU() const
int waferV() const
int waferV() const
Definition: HFNoseDetId.h:78
int32_t waferV(const int32_t index)

◆ isHalfCell()

bool RecHitTools::isHalfCell ( const DetId id) const

Definition at line 415 of file RecHitTools.cc.

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

415  {
416  bool ishalf = false;
417  if (id.det() == DetId::Forward) {
418  HGCalDetId hid(id);
419  auto geom = getSubdetectorGeometry(hid);
420  auto ddd = get_ddd(geom, hid);
421  const int waferType = ddd->waferTypeT(hid.waferType());
422  return ddd->isHalfCell(waferType, hid.cell());
423  }
424  //new geometry is always false
425  return ishalf;
426 }
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
Definition: RecHitTools.cc:119

◆ isOnlySilicon()

bool RecHitTools::isOnlySilicon ( const unsigned int  layer) const

Definition at line 435 of file RecHitTools.cc.

References bhLastLayer_, bhOffset_, and phase1PixelTopology::layer.

435  {
436  // HFnose TODO
437  bool isonlysilicon = (layer % bhLastLayer_) < bhOffset_;
438  return isonlysilicon;
439 }
unsigned int bhLastLayer_
Definition: RecHitTools.h:90
unsigned int bhOffset_
Definition: RecHitTools.h:90
constexpr std::array< uint8_t, layerIndexSize > layer

◆ isScintillator()

bool RecHitTools::isScintillator ( const DetId id) const

Definition at line 433 of file RecHitTools.cc.

References DetId::HGCalHSc.

Referenced by getScintDEtaDPhi().

433 { return id.det() == DetId::HGCalHSc; }

◆ isSilicon()

bool RecHitTools::isSilicon ( const DetId id) const

Definition at line 428 of file RecHitTools.cc.

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

Referenced by FWRecoGeometryESProducer::addCaloGeometry(), ticl::ClusterFilterByAlgoAndSize::filter(), and ticl::ClusterFilterByAlgoAndSizeAndLayerRange::filter().

428  {
429  return (id.det() == DetId::HGCalEE || id.det() == DetId::HGCalHSi ||
430  (id.det() == DetId::Forward && id.subdetId() == static_cast<int>(HFNose)));
431 }

◆ lastLayer()

unsigned int hgcal::RecHitTools::lastLayer ( bool  nose = false) const
inline

Definition at line 79 of file RecHitTools.h.

References bhLastLayer_, and noseLastLayer_.

Referenced by HGCalClusteringAlgoBase::getEventSetup(), TICLLayerTileProducer::produce(), and setGeometry().

79 { return (nose ? noseLastLayer_ : bhLastLayer_); }
unsigned int bhLastLayer_
Definition: RecHitTools.h:90
unsigned int noseLastLayer_
Definition: RecHitTools.h:90

◆ lastLayerBH()

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

Definition at line 78 of file RecHitTools.h.

References bhLastLayer_.

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

78 { return bhLastLayer_; }
unsigned int bhLastLayer_
Definition: RecHitTools.h:90

◆ lastLayerEE()

unsigned int hgcal::RecHitTools::lastLayerEE ( bool  nose = false) const
inline

◆ lastLayerFH()

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

Definition at line 76 of file RecHitTools.h.

References fhLastLayer_.

Referenced by RealisticSimClusterMapper::buildClusters(), and HGCalClusteringAlgoBase::getEventSetup().

76 { return fhLastLayer_; }
unsigned int fhLastLayer_
Definition: RecHitTools.h:90

◆ maskCell()

bool RecHitTools::maskCell ( const DetId id,
int  corners = 3 
) const

Definition at line 488 of file RecHitTools.cc.

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

488  {
489  if (id.det() == DetId::Hcal) {
490  return false;
491  } else {
492  auto hg = static_cast<const HGCalGeometry*>(getSubdetectorGeometry(id));
493  return hg->topology().dddConstants().maskCell(id, corners);
494  }
495 }
bool maskCell(const DetId &id, int corners) const
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
Definition: RecHitTools.cc:119
const HGCalTopology & topology() const
const HGCalDDDConstants & dddConstants() const
Definition: HGCalTopology.h:98

◆ maxNumberOfWafersPerLayer()

unsigned int hgcal::RecHitTools::maxNumberOfWafersPerLayer ( bool  nose = false) const
inline

Definition at line 81 of file RecHitTools.h.

References maxNumberOfWafersNose_, and maxNumberOfWafersPerLayer_.

81  {
83  }
unsigned int maxNumberOfWafersNose_
Definition: RecHitTools.h:91
unsigned int maxNumberOfWafersPerLayer_
Definition: RecHitTools.h:91

◆ setGeometry()

void RecHitTools::setGeometry ( CaloGeometry const &  geom)

Definition at line 68 of file RecHitTools.cc.

References bhFirstLayer_, bhLastLayer_, bhMaxIphi_, bhOffset_, eeOffset_, fhLastLayer_, fhOffset_, DetId::Forward, ForwardEmpty, relativeConstraints::geom, geom_, geometryType_, CaloGeometry::getSubdetectorGeometry(), DetId::Hcal, HcalEndcap, HFNose, DetId::HGCalEE, DetId::HGCalHSc, DetId::HGCalHSi, HGCEE, HGCHEF, lastLayer(), hgcalTopologyTester_cfi::layers, SiStripPI::max, maxNumberOfWafersNose_, maxNumberOfWafersPerLayer_, and noseLastLayer_.

Referenced by FWRecoGeometryESProducer::addCaloGeometry(), HeterogeneousHGCalRecHitsValidator::analyze(), HGCalShowerSeparation::analyze(), HGCalHitCalibration::analyze(), PositionAtECalEntranceComputer::beginEvent(), TICLLayerTileProducer::beginRun(), FilteredLayerClustersProducer::beginRun(), TICLTrackstersEdgesValidation::dqmBeginRun(), HGCalEgammaIDHelper::eventInit(), hgcal::ClusterTools::getEventSetup(), HGCalDepthPreClusterer::getEventSetup(), HGCal3DClustering::getEventSetup(), HGCalClusteringAlgoBase::getEventSetup(), PFHGCalRecHitCreator< DET, Layer, det, subdet >::importRecHits(), TrackstersMergeProducer::produce(), SimTrackstersProducer::produce(), HGCalRecHitAbsAlgo::set(), and RealisticSimClusterMapper::update().

68  {
69  geom_ = &geom;
70  unsigned int wmaxEE(0), wmaxFH(0);
71  auto geomEE = static_cast<const HGCalGeometry*>(
73  //check if it's the new geometry
74  if (geomEE) {
75  geometryType_ = 1;
76  eeOffset_ = (geomEE->topology().dddConstants()).getLayerOffset();
77  wmaxEE = (geomEE->topology().dddConstants()).waferCount(0);
78  auto geomFH = static_cast<const HGCalGeometry*>(
80  fhOffset_ = (geomFH->topology().dddConstants()).getLayerOffset();
81  wmaxFH = (geomFH->topology().dddConstants()).waferCount(0);
82  fhLastLayer_ = fhOffset_ + (geomFH->topology().dddConstants()).lastLayer(true);
83  auto geomBH = static_cast<const HGCalGeometry*>(
85  bhOffset_ = (geomBH->topology().dddConstants()).getLayerOffset();
86  bhFirstLayer_ = bhOffset_ + (geomBH->topology().dddConstants()).firstLayer();
87  bhLastLayer_ = bhOffset_ + (geomBH->topology().dddConstants()).lastLayer(true);
88  bhMaxIphi_ = geomBH->topology().dddConstants().maxCells(true);
89  } else {
90  geometryType_ = 0;
91  geomEE =
93  eeOffset_ = (geomEE->topology().dddConstants()).getLayerOffset();
94  wmaxEE = 1 + (geomEE->topology().dddConstants()).waferMax();
95  auto geomFH =
97  fhOffset_ = (geomFH->topology().dddConstants()).getLayerOffset();
98  fhLastLayer_ = fhOffset_ + (geomFH->topology().dddConstants()).layers(true);
100  bhFirstLayer_ = bhOffset_ + 1;
101  wmaxFH = 1 + (geomFH->topology().dddConstants()).waferMax();
102  auto geomBH =
104  bhLastLayer_ = bhOffset_ + (geomBH->topology().dddConstants())->getMaxDepth(1);
105  }
106  maxNumberOfWafersPerLayer_ = std::max(wmaxEE, wmaxFH);
107  // For nose geometry
108  auto geomNose =
110  if (geomNose) {
111  maxNumberOfWafersNose_ = (geomNose->topology().dddConstants()).waferCount(0);
112  noseLastLayer_ = (geomNose->topology().dddConstants()).layers(true);
113  } else {
115  noseLastLayer_ = 0;
116  }
117 }
unsigned int bhFirstLayer_
Definition: RecHitTools.h:90
unsigned int fhLastLayer_
Definition: RecHitTools.h:90
unsigned int bhLastLayer_
Definition: RecHitTools.h:90
unsigned int maxNumberOfWafersNose_
Definition: RecHitTools.h:91
unsigned int fhOffset_
Definition: RecHitTools.h:90
unsigned int bhOffset_
Definition: RecHitTools.h:90
unsigned int maxNumberOfWafersPerLayer_
Definition: RecHitTools.h:91
const CaloGeometry * geom_
Definition: RecHitTools.h:89
unsigned int eeOffset_
Definition: RecHitTools.h:90
unsigned int noseLastLayer_
Definition: RecHitTools.h:90
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:34
unsigned int lastLayer(bool nose=false) const
Definition: RecHitTools.h:79

◆ zside()

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 FWRecoGeometryESProducer::addCaloGeometry(), and 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) {
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 }
constexpr int zside() const
get the z-side of the cell (1/-1)
Definition: HcalDetId.h:141
int zside(const DetId &id) const
Definition: RecHitTools.cc:163
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)
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: HFNoseDetId.h:53

Member Data Documentation

◆ bhFirstLayer_

unsigned int hgcal::RecHitTools::bhFirstLayer_
private

Definition at line 90 of file RecHitTools.h.

Referenced by firstAndLastLayer(), firstLayerBH(), and setGeometry().

◆ bhLastLayer_

unsigned int hgcal::RecHitTools::bhLastLayer_
private

Definition at line 90 of file RecHitTools.h.

Referenced by firstAndLastLayer(), isOnlySilicon(), lastLayer(), lastLayerBH(), and setGeometry().

◆ bhMaxIphi_

int hgcal::RecHitTools::bhMaxIphi_
private

Definition at line 93 of file RecHitTools.h.

Referenced by getScintMaxIphi(), and setGeometry().

◆ bhOffset_

unsigned int hgcal::RecHitTools::bhOffset_
private

Definition at line 90 of file RecHitTools.h.

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

◆ eeOffset_

unsigned int hgcal::RecHitTools::eeOffset_
private

Definition at line 90 of file RecHitTools.h.

Referenced by firstAndLastLayer(), and setGeometry().

◆ fhLastLayer_

unsigned int hgcal::RecHitTools::fhLastLayer_
private

Definition at line 90 of file RecHitTools.h.

Referenced by firstAndLastLayer(), lastLayerFH(), and setGeometry().

◆ fhOffset_

unsigned int hgcal::RecHitTools::fhOffset_
private

◆ geom_

const CaloGeometry* hgcal::RecHitTools::geom_
private

◆ geometryType_

int hgcal::RecHitTools::geometryType_
private

Definition at line 92 of file RecHitTools.h.

Referenced by getGeometryType(), getPositionLayer(), and setGeometry().

◆ maxNumberOfWafersNose_

unsigned int hgcal::RecHitTools::maxNumberOfWafersNose_
private

Definition at line 91 of file RecHitTools.h.

Referenced by maxNumberOfWafersPerLayer(), and setGeometry().

◆ maxNumberOfWafersPerLayer_

unsigned int hgcal::RecHitTools::maxNumberOfWafersPerLayer_
private

Definition at line 91 of file RecHitTools.h.

Referenced by maxNumberOfWafersPerLayer(), and setGeometry().

◆ noseLastLayer_

unsigned int hgcal::RecHitTools::noseLastLayer_
private

Definition at line 90 of file RecHitTools.h.

Referenced by firstAndLastLayer(), lastLayer(), and setGeometry().