CMS 3D CMS Logo

List of all members | Classes | Public Member Functions | Static Public Attributes | Private Member Functions | Private Attributes
HGCalTopology Class Reference

#include <HGCalTopology.h>

Inheritance diagram for HGCalTopology:
CaloSubdetectorTopology

Classes

struct  DecodedDetId
 

Public Member Functions

unsigned int allGeomModules () const
 
const HGCalDDDConstantsdddConstants () const
 
DecodedDetId decode (const DetId &id) const
 
DetId denseId2detId (uint32_t denseId) const override
 
DetId::Detector detector () const
 
bool detectorType () const
 
virtual uint32_t detId2denseGeomId (const DetId &id) const
 
uint32_t detId2denseId (const DetId &id) const override
 Dense indexing. More...
 
std::vector< DetIddown (const DetId &id) const override
 
std::vector< DetIdeast (const DetId &id) const override
 
DetId encode (const DecodedDetId &id_) const
 
DecodedDetId geomDenseId2decId (const uint32_t &hi) const
 
HGCalGeometryMode::GeometryMode geomMode () const
 Geometry mode. More...
 
DetId goEast (const DetId &id) const override
 move the Topology east (positive ix) More...
 
DetId goNorth (const DetId &id) const override
 move the Topology north (increment iy) More...
 
DetId goSouth (const DetId &id) const override
 move the Topology south (decrement iy) More...
 
DetId goWest (const DetId &id) const override
 move the Topology west (negative ix) More...
 
 HGCalTopology (const HGCalDDDConstants &hdcons, int subdet)
 create a new Topology More...
 
bool isHFNose () const
 
bool maskCell (const DetId &id, int corners=3) const
 
std::vector< DetIdneighbors (const DetId &id) const
 
std::vector< DetIdnorth (const DetId &id) const override
 
DetId offsetBy (const DetId startId, int nrStepsX, int nrStepsY) const
 
std::vector< DetIdsouth (const DetId &id) const override
 
ForwardSubdetector subDetector () const
 
DetId switchZSide (const DetId startId) const
 
bool tileTrapezoid () const
 
unsigned int totalGeomModules () const
 
unsigned int totalModules () const
 
std::vector< DetIdup (const DetId &id) const override
 
bool valid (const DetId &id) const override
 Is this a valid cell id. More...
 
bool valid (const DetId &id, int cornerMin) const
 
bool validHashIndex (uint32_t ix) const
 
bool validModule (const DetId &id, int cornerMin) const
 
bool waferHexagon6 () const
 
bool waferHexagon8 () const
 
std::vector< DetIdwest (const DetId &id) const override
 
 ~HGCalTopology () override=default
 virtual destructor More...
 
- Public Member Functions inherited from CaloSubdetectorTopology
 CaloSubdetectorTopology ()
 standard constructor More...
 
virtual DetId denseId2detId (unsigned int) const
 return a linear packed id More...
 
virtual bool denseIdConsistent (int topoVer) const
 return whether this topology is consistent with the numbering in the given topology More...
 
virtual std::vector< DetIdgetAllNeighbours (const DetId &id) const
 
virtual std::vector< DetIdgetNeighbours (const DetId &id, const CaloDirection &dir) const
 
virtual std::vector< DetIdgetWindow (const DetId &id, const int &northSouthSize, const int &eastWestSize) const
 
virtual DetId goDown (const DetId &id) const
 
virtual DetId goUp (const DetId &id) const
 
virtual unsigned int ncells () const
 return a count of valid cells (for dense indexing use) More...
 
virtual int topoVersion () const
 return a version which identifies the given topology More...
 
virtual ~CaloSubdetectorTopology ()
 virtual destructor More...
 

Static Public Attributes

static const int subSectors_ = 2
 Use subSector in square mode as wafer type in hexagon mode. More...
 

Private Member Functions

void addHGCSCintillatorId (std::vector< DetId > &ids, int zside, int type, int lay, int iradius, int iphi) const
 add DetId of Scintillator and Silicon type if valid More...
 
void addHGCSiliconId (std::vector< DetId > &ids, int det, int zside, int type, int lay, int waferU, int waferV, int cellU, int cellV) const
 
DetId changeXY (const DetId &id, int nrStepsX, int nrStepsY) const
 move the nagivator along x, y More...
 
DetId changeZ (const DetId &id, int nrStepsZ) const
 move the nagivator along z More...
 

Private Attributes

int cellMax_
 
int cells_
 
DetId::Detector det_
 
int firstLay_
 
const HGCalDDDConstantshdcons_
 
int kHGeomHalf_
 
int kHGhalf_
 
int kHGhalfType_
 
unsigned int kSizeForDenseIndexing
 
int layers_
 
HGCalGeometryMode::GeometryMode mode_
 
int sectors_
 
ForwardSubdetector subdet_
 
int types_
 
int waferMax_
 
int waferOff_
 

Additional Inherited Members

- Protected Types inherited from CaloSubdetectorTopology
typedef std::pair< int, int > Coordinate
 
- Protected Member Functions inherited from CaloSubdetectorTopology
Coordinate getNeighbourIndex (const Coordinate &coord, const CaloDirection &dir) const
 

Detailed Description

Definition at line 12 of file HGCalTopology.h.

Constructor & Destructor Documentation

◆ HGCalTopology()

HGCalTopology::HGCalTopology ( const HGCalDDDConstants hdcons,
int  subdet 
)

create a new Topology

Definition at line 13 of file HGCalTopology.cc.

References cellMax_, cells_, det_, firstLay_, HGCalDDDConstants::firstLayer(), DetId::Forward, ForwardEmpty, HGCalDDDConstants::geomMode(), hdcons_, HFNose, kHGeomHalf_, kHGhalf_, kHGhalfType_, kSizeForDenseIndexing, HGCalDDDConstants::layers(), layers_, HGCalDDDConstants::maxCells(), HGCalDDDConstants::maxCellUV(), mode_, HGCalDDDConstants::sectors(), sectors_, subdet_, tileTrapezoid(), types_, waferHexagon6(), waferMax_, waferOff_, and HGCalDDDConstants::waferUVMax().

13  : hdcons_(hdcons) {
15  layers_ = hdcons_.layers(true);
16  cells_ = hdcons_.maxCells(true);
20  waferMax_ = 2 * waferOff_ + 1;
23  if (waferHexagon6()) {
25  subdet_ = (ForwardSubdetector)(det);
27  types_ = 2;
28  } else if (det == (int)(DetId::Forward)) {
30  subdet_ = HFNose;
32  types_ = 3;
33  } else if (tileTrapezoid()) {
34  det_ = (DetId::Detector)(det);
37  types_ = 3;
38  } else {
39  det_ = (DetId::Detector)(det);
42  types_ = 3;
43  }
45  kSizeForDenseIndexing = static_cast<unsigned int>(2 * kHGhalf_);
46 #ifdef EDM_ML_DEBUG
47  edm::LogVerbatim("HGCalGeom") << "HGCalTopology initialized for detector " << det << ":" << det_ << ":" << subdet_
48  << " having " << sectors_ << " Sectors, " << layers_ << " Layers from " << firstLay_
49  << ", " << cells_ << " cells and total channels " << kSizeForDenseIndexing << ":"
50  << (2 * kHGeomHalf_);
51 #endif
52 }
Log< level::Info, true > LogVerbatim
HGCalGeometryMode::GeometryMode mode_
bool tileTrapezoid() const
bool waferHexagon6() const
int firstLayer() const
HGCalGeometryMode::GeometryMode geomMode() const
ForwardSubdetector
DetId::Detector det_
int maxCells(bool reco) const
unsigned int layers(bool reco) const
ForwardSubdetector subdet_
const HGCalDDDConstants & hdcons_
Detector
Definition: DetId.h:24
unsigned int kSizeForDenseIndexing

◆ ~HGCalTopology()

HGCalTopology::~HGCalTopology ( )
overridedefault

virtual destructor

Member Function Documentation

◆ addHGCSCintillatorId()

void HGCalTopology::addHGCSCintillatorId ( std::vector< DetId > &  ids,
int  zside,
int  type,
int  lay,
int  iradius,
int  iphi 
) const
private

add DetId of Scintillator and Silicon type if valid

Definition at line 565 of file HGCalTopology.cc.

References hdcons_, l1ctLayer2EG_cff::id, LEDCalibrationChannels::iphi, HGCalDDDConstants::isValidTrap(), and ecaldqm::zside().

Referenced by neighbors().

566  {
567 #ifdef EDM_ML_DEBUG
568  edm::LogVerbatim("HGCalGeom") << "addHGCSCintillatorId " << zside << ":" << type << ":" << lay << ":" << iradius
569  << ":" << iphi << " ==> Validity " << hdcons_.isValidTrap(zside, lay, iradius, iphi);
570 #endif
571  if (hdcons_.isValidTrap(zside, lay, iradius, iphi)) {
572  HGCScintillatorDetId id(type, lay, zside * iradius, iphi);
573  ids.emplace_back(DetId(id));
574  }
575 }
Log< level::Info, true > LogVerbatim
int zside(DetId const &)
bool isValidTrap(int zside, int lay, int ieta, int iphi) const
const HGCalDDDConstants & hdcons_
Definition: DetId.h:17

◆ addHGCSiliconId()

void HGCalTopology::addHGCSiliconId ( std::vector< DetId > &  ids,
int  det,
int  zside,
int  type,
int  lay,
int  waferU,
int  waferV,
int  cellU,
int  cellV 
) const
private

Definition at line 577 of file HGCalTopology.cc.

References hdcons_, l1ctLayer2EG_cff::id, HGCalDDDConstants::isValidHex8(), HGCalWaferIndex::waferU(), HGCalWaferIndex::waferV(), and ecaldqm::zside().

Referenced by neighbors().

578  {
579 #ifdef EDM_ML_DEBUG
580  edm::LogVerbatim("HGCalGeom") << "addHGCSiliconId " << det << ":" << zside << ":" << type << ":" << lay << ":"
581  << waferU << ":" << waferV << ":" << cellU << ":" << cellV << " ==> Validity "
582  << hdcons_.isValidHex8(lay, waferU, waferV, cellU, cellV);
583 #endif
584  if (hdcons_.isValidHex8(lay, waferU, waferV, cellU, cellV)) {
585  HGCSiliconDetId id((DetId::Detector)(det), zside, type, lay, waferU, waferV, cellU, cellV);
586  ids.emplace_back(DetId(id));
587  }
588 }
Log< level::Info, true > LogVerbatim
int32_t waferU(const int32_t index)
bool isValidHex8(int lay, int waferU, int waferV, bool fullAndPart=false) const
int zside(DetId const &)
const HGCalDDDConstants & hdcons_
Definition: DetId.h:17
Detector
Definition: DetId.h:24
int32_t waferV(const int32_t index)

◆ allGeomModules()

unsigned int HGCalTopology::allGeomModules ( ) const

Definition at line 54 of file HGCalTopology.cc.

References hdcons_, HGCalDDDConstants::numberCells(), tileTrapezoid(), and HGCalDDDConstants::wafers().

Referenced by HGCalGeometryLoader::build().

54  {
55  return (tileTrapezoid() ? (unsigned int)(2 * hdcons_.numberCells(true)) : (unsigned int)(2 * hdcons_.wafers()));
56 }
bool tileTrapezoid() const
const HGCalDDDConstants & hdcons_
int numberCells(bool reco) const

◆ changeXY()

DetId HGCalTopology::changeXY ( const DetId id,
int  nrStepsX,
int  nrStepsY 
) const
private

move the nagivator along x, y

Definition at line 664 of file HGCalTopology.cc.

Referenced by goEast(), goNorth(), goSouth(), goWest(), and offsetBy().

664 { return DetId(); }
Definition: DetId.h:17

◆ changeZ()

DetId HGCalTopology::changeZ ( const DetId id,
int  nrStepsZ 
) const
private

move the nagivator along z

Definition at line 666 of file HGCalTopology.cc.

Referenced by down(), and up().

666 { return DetId(); }
Definition: DetId.h:17

◆ dddConstants()

const HGCalDDDConstants& HGCalTopology::dddConstants ( ) const
inline

◆ decode()

HGCalTopology::DecodedDetId HGCalTopology::decode ( const DetId id) const

Definition at line 590 of file HGCalTopology.cc.

References det_, DetId::Forward, HFNose, l1ctLayer2EG_cff::id, heavyIonCSV_trainingSettings::idx, createfilelist::int, subdet_, tileTrapezoid(), and waferHexagon6().

Referenced by detId2denseGeomId(), detId2denseId(), HGCalGeometry::get8Corners(), HGCalGeometry::getClosestCell(), HGCalGeometry::getClosestCellHex(), HGCalGeometry::getCorners(), HGCalGeometry::getNewCorners(), HGCalGeometry::getPosition(), neighbors(), HGCalGeometry::neighborZ(), HGCalGeometry::newCell(), switchZSide(), valid(), and validModule().

590  {
592  if (waferHexagon6()) {
593  HGCalDetId id(startId);
594  idx.iCell1 = id.cell();
595  idx.iCell2 = 0;
596  idx.iLay = id.layer();
597  idx.iSec1 = id.wafer();
598  idx.iSec2 = 0;
599  idx.iType = id.waferType();
600  idx.zSide = id.zside();
601  idx.det = id.subdetId();
602  } else if (tileTrapezoid()) {
603  HGCScintillatorDetId id(startId);
604  idx.iCell1 = id.iphi();
605  idx.iCell2 = 0;
606  idx.iLay = id.layer();
607  idx.iSec1 = id.ietaAbs();
608  idx.iSec2 = 0;
609  idx.iType = id.type();
610  idx.zSide = id.zside();
611  idx.det = (int)(id.subdet());
613  HFNoseDetId id(startId);
614  idx.iCell1 = id.cellU();
615  idx.iCell2 = id.cellV();
616  idx.iLay = id.layer();
617  idx.iSec1 = id.waferU();
618  idx.iSec2 = id.waferV();
619  idx.iType = id.type();
620  idx.zSide = id.zside();
621  idx.det = (int)(id.subdet());
622  } else {
623  HGCSiliconDetId id(startId);
624  idx.iCell1 = id.cellU();
625  idx.iCell2 = id.cellV();
626  idx.iLay = id.layer();
627  idx.iSec1 = id.waferU();
628  idx.iSec2 = id.waferV();
629  idx.iType = id.type();
630  idx.zSide = id.zside();
631  idx.det = (int)(id.subdet());
632  }
633  return idx;
634 }
bool tileTrapezoid() const
bool waferHexagon6() const
DetId::Detector det_
ForwardSubdetector subdet_

◆ denseId2detId()

DetId HGCalTopology::denseId2detId ( uint32_t  denseId) const
override

Definition at line 389 of file HGCalTopology.cc.

References cellMax_, encode(), firstLay_, l1ctLayer2EG_cff::id, createfilelist::int, kHGhalfType_, layers_, sectors_, tileTrapezoid(), types_, validHashIndex(), waferHexagon6(), waferMax_, and waferOff_.

389  {
391  if (validHashIndex(hi)) {
392  id.zSide = ((int)(hi) < kHGhalfType_ ? -1 : 1);
393  int di = ((int)(hi) % kHGhalfType_);
394  if (waferHexagon6()) {
395  int type = (di % types_);
396  id.iType = (type == 0 ? -1 : 1);
397  id.iSec1 = (((di - type) / types_) % sectors_);
398  id.iLay = (((((di - type) / types_) - id.iSec1 + 1) / sectors_) % layers_ + 1);
399  id.iCell1 = (((((di - type) / types_) - id.iSec1 + 1) / sectors_ - id.iLay + 1) / layers_ + 1);
400 #ifdef EDM_ML_DEBUG
401  edm::LogVerbatim("HGCalGeom") << "Input Hex " << hi << " o/p " << id.zSide << ":" << id.iLay << ":" << id.iType
402  << ":" << id.iSec1 << ":" << id.iCell1;
403 #endif
404  } else if (tileTrapezoid()) {
405  int type = (di % types_);
406  id.iType = type;
407  id.iSec1 = (((di - type) / types_) % sectors_) + 1;
408  id.iLay = (((((di - type) / types_) - id.iSec1 + 1) / sectors_) % layers_ + firstLay_);
409  id.iCell1 = (((((di - type) / types_) - id.iSec1 + 1) / sectors_ - id.iLay + firstLay_) / layers_ + 1);
410 #ifdef EDM_ML_DEBUG
411  edm::LogVerbatim("HGCalGeom") << "Input Trap " << hi << " o/p " << id.zSide << ":" << id.iLay << ":" << id.iType
412  << ":" << id.iSec1 << ":" << id.iCell1;
413 #endif
414  } else {
415  int type = (di % types_);
416  id.iType = type;
417  di = (di - type) / types_;
418  id.iSec2 = (di % waferMax_) - waferOff_;
419  di = (di - id.iSec2 - waferOff_) / waferMax_;
420  id.iSec1 = (di % waferMax_) - waferOff_;
421  di = (di - id.iSec1 - waferOff_) / waferMax_;
422  id.iLay = (di % layers_) + 1;
423  di = (di - id.iLay + 1) / layers_;
424  id.iCell2 = (di % cellMax_);
425  id.iCell1 = (di - id.iCell2) / cellMax_;
426 #ifdef EDM_ML_DEBUG
427  edm::LogVerbatim("HGCalGeom") << "Input Hex8 " << hi << " o/p " << id.zSide << ":" << id.iLay << ":" << id.iType
428  << ":" << id.iSec1 << ":" << id.iSec2 << ":" << id.iCell1 << ":" << id.iCell2;
429 #endif
430  }
431  }
432  return encode(id);
433 }
Log< level::Info, true > LogVerbatim
bool tileTrapezoid() const
bool waferHexagon6() const
Definition: EPCuts.h:4
bool validHashIndex(uint32_t ix) const
Definition: HGCalTopology.h:89
DetId encode(const DecodedDetId &id_) const

◆ detector()

DetId::Detector HGCalTopology::detector ( ) const
inline

Definition at line 117 of file HGCalTopology.h.

References det_.

Referenced by HGCalGeometryLoader::build(), HGCalGeometry::getClosestCell(), and HGCalGeometry::getClosestCellHex().

117 { return det_; }
DetId::Detector det_

◆ detectorType()

bool HGCalTopology::detectorType ( ) const
inline

Definition at line 119 of file HGCalTopology.h.

119 { return false; }

◆ detId2denseGeomId()

uint32_t HGCalTopology::detId2denseGeomId ( const DetId id) const
virtual

Definition at line 435 of file HGCalTopology.cc.

References cellMax_, decode(), firstLay_, heavyIonCSV_trainingSettings::idx, kHGeomHalf_, layers_, sectors_, tileTrapezoid(), waferHexagon6(), waferMax_, and waferOff_.

Referenced by HGCalGeometry::getGeometry(), HGCalGeometry::getSummary(), HGCalGeometry::indexFor(), HGCalGeometry::newCell(), and HGCalGeometry::present().

435  {
437  uint32_t idx;
438  if (waferHexagon6()) {
439  idx = (uint32_t)(((id.zSide > 0) ? kHGeomHalf_ : 0) + (id.iLay - 1) * sectors_ + id.iSec1);
440 #ifdef EDM_ML_DEBUG
441  edm::LogVerbatim("HGCalGeom") << "Geom Hex I/P " << id.zSide << ":" << id.iLay << ":" << id.iSec1 << ":" << id.iType
442  << " Constants " << kHGeomHalf_ << ":" << layers_ << ":" << sectors_ << " o/p "
443  << idx;
444 #endif
445  } else if (tileTrapezoid()) {
446  idx = (uint32_t)(((id.zSide > 0) ? kHGeomHalf_ : 0) +
447  (((id.iLay - firstLay_) * sectors_ + id.iSec1 - 1) * cellMax_ + id.iCell1 - 1));
448 #ifdef EDM_ML_DEBUG
449  edm::LogVerbatim("HGCalGeom") << "Geom Trap I/P " << id.zSide << ":" << id.iLay << ":" << id.iSec1 << ":"
450  << id.iCell1 << ":" << id.iType << " Constants " << kHGeomHalf_ << ":" << layers_
451  << ":" << firstLay_ << ":" << sectors_ << ":" << cellMax_ << " o/p " << idx;
452 #endif
453  } else {
454  idx = (uint32_t)(((id.zSide > 0) ? kHGeomHalf_ : 0) +
455  (((id.iLay - 1) * waferMax_ + id.iSec1 + waferOff_) * waferMax_ + id.iSec2 + waferOff_));
456 #ifdef EDM_ML_DEBUG
457  edm::LogVerbatim("HGCalGeom") << "Geom Hex8 I/P " << id.zSide << ":" << id.iLay << ":" << id.iSec1 << ":"
458  << id.iSec2 << ":" << id.iType << " Constants " << kHGeomHalf_ << ":" << layers_
459  << ":" << waferMax_ << ":" << waferOff_ << " o/p " << idx;
460 #endif
461  }
462  return idx;
463 }
Log< level::Info, true > LogVerbatim
bool tileTrapezoid() const
bool waferHexagon6() const
DecodedDetId decode(const DetId &id) const

◆ detId2denseId()

uint32_t HGCalTopology::detId2denseId ( const DetId id) const
overridevirtual

Dense indexing.

Reimplemented from CaloSubdetectorTopology.

Definition at line 350 of file HGCalTopology.cc.

References cellMax_, decode(), firstLay_, heavyIonCSV_trainingSettings::idx, kHGeomHalf_, kHGhalfType_, layers_, sectors_, tileTrapezoid(), types_, waferHexagon6(), waferMax_, and waferOff_.

350  {
352  uint32_t idx;
353  if (waferHexagon6()) {
354  int type = (id.iType > 0) ? 1 : 0;
355  idx = (uint32_t)(((id.zSide > 0) ? kHGhalfType_ : 0) +
356  ((((id.iCell1 - 1) * layers_ + id.iLay - 1) * sectors_ + id.iSec1) * types_ + type));
357 #ifdef EDM_ML_DEBUG
358  edm::LogVerbatim("HGCalGeom") << "Input Hex " << id.zSide << ":" << id.iLay << ":" << id.iSec1 << ":" << id.iCell1
359  << ":" << id.iType << " Constants " << kHGeomHalf_ << ":" << layers_ << ":"
360  << sectors_ << ":" << types_ << " o/p " << idx;
361 #endif
362  } else if (tileTrapezoid()) {
363  idx =
364  (uint32_t)(((id.zSide > 0) ? kHGhalfType_ : 0) +
365  ((((id.iCell1 - 1) * layers_ + id.iLay - firstLay_) * sectors_ + id.iSec1 - 1) * types_ + id.iType));
366 #ifdef EDM_ML_DEBUG
367  edm::LogVerbatim("HGCalGeom") << "Input Trap " << id.zSide << ":" << id.iLay << ":" << id.iSec1 << ":" << id.iCell1
368  << ":" << id.iType << " Constants " << kHGeomHalf_ << ":" << layers_ << ":"
369  << sectors_ << ":" << types_ << " o/p " << idx;
370 #endif
371  } else {
372  idx =
373  (uint32_t)(((id.zSide > 0) ? kHGhalfType_ : 0) +
374  (((((id.iCell1 * cellMax_ + id.iCell2) * layers_ + id.iLay - 1) * waferMax_ + id.iSec1 + waferOff_) *
375  waferMax_ +
376  id.iSec2 + waferOff_) *
377  types_ +
378  id.iType));
379 #ifdef EDM_ML_DEBUG
380  edm::LogVerbatim("HGCalGeom") << "Input Hex8 " << id.zSide << ":" << id.iLay << ":" << id.iSec1 << ":" << id.iSec2
381  << ":" << id.iCell1 << ":" << id.iCell2 << ":" << id.iType << " Constants "
382  << kHGeomHalf_ << ":" << cellMax_ << ":" << layers_ << ":" << waferMax_ << ":"
383  << waferOff_ << ":" << types_ << " o/p " << idx;
384 #endif
385  }
386  return idx;
387 }
Log< level::Info, true > LogVerbatim
bool tileTrapezoid() const
bool waferHexagon6() const
DecodedDetId decode(const DetId &id) const

◆ down()

std::vector<DetId> HGCalTopology::down ( const DetId id) const
inlineoverridevirtual

Get the neighbors of the given cell in down direction (inward)

Implements CaloSubdetectorTopology.

Definition at line 68 of file HGCalTopology.h.

References changeZ().

68  {
69  DetId nextId = changeZ(id, -1);
70  std::vector<DetId> vNeighborsDetId;
71  if (!(nextId == DetId(0)))
72  vNeighborsDetId.emplace_back(nextId);
73  return vNeighborsDetId;
74  }
Definition: DetId.h:17
DetId changeZ(const DetId &id, int nrStepsZ) const
move the nagivator along z

◆ east()

std::vector<DetId> HGCalTopology::east ( const DetId id) const
inlineoverridevirtual

Get the neighbors of the given cell in east direction

Implements CaloSubdetectorTopology.

Definition at line 42 of file HGCalTopology.h.

References goEast().

42  {
43  DetId nextId = goEast(id);
44  std::vector<DetId> vNeighborsDetId;
45  if (!(nextId == DetId(0)))
46  vNeighborsDetId.emplace_back(nextId);
47  return vNeighborsDetId;
48  }
DetId goEast(const DetId &id) const override
move the Topology east (positive ix)
Definition: HGCalTopology.h:41
Definition: DetId.h:17

◆ encode()

DetId HGCalTopology::encode ( const DecodedDetId id_) const

Definition at line 636 of file HGCalTopology.cc.

References det_, DetId::Forward, hdcons_, HFNose, l1ctLayer2EG_cff::id, heavyIonCSV_trainingSettings::idx, HGCScintillatorDetId::layer(), DetId::rawId(), HGCScintillatorDetId::ring(), HGCScintillatorDetId::setSiPM(), HGCScintillatorDetId::setType(), subdet_, tileTrapezoid(), HGCalDDDConstants::tileType(), and waferHexagon6().

Referenced by denseId2detId(), HGCalGeometry::getClosestCell(), HGCalGeometry::getClosestCellHex(), HGCalGeometry::newCell(), and switchZSide().

636  {
637  DetId id;
638 #ifdef EDM_ML_DEBUG
639  edm::LogVerbatim("HGCalGeomX") << "Encode " << idx.det << ":" << idx.zSide << ":" << idx.iType << ":" << idx.iLay
640  << ":" << idx.iSec1 << ":" << idx.iSec2 << ":" << idx.iCell1 << ":" << idx.iCell2;
641 #endif
642  if (waferHexagon6()) {
643  id =
644  HGCalDetId((ForwardSubdetector)(idx.det), idx.zSide, idx.iLay, ((idx.iType > 0) ? 1 : 0), idx.iSec1, idx.iCell1)
645  .rawId();
646  } else if (tileTrapezoid()) {
647  HGCScintillatorDetId hid(idx.iType, idx.iLay, idx.zSide * idx.iSec1, idx.iCell1);
648  std::pair<int, int> typm = hdcons_.tileType(hid.layer(), hid.ring(), 0);
649  if (typm.first >= 0) {
650  hid.setType(typm.first);
651  hid.setSiPM(typm.second);
652  }
653  id = hid.rawId();
655  id = HFNoseDetId(idx.zSide, idx.iType, idx.iLay, idx.iSec1, idx.iSec2, idx.iCell1, idx.iCell2).rawId();
656  } else {
657  id = HGCSiliconDetId(
658  (DetId::Detector)(idx.det), idx.zSide, idx.iType, idx.iLay, idx.iSec1, idx.iSec2, idx.iCell1, idx.iCell2)
659  .rawId();
660  }
661  return id;
662 }
Log< level::Info, true > LogVerbatim
bool tileTrapezoid() const
bool waferHexagon6() const
ForwardSubdetector
DetId::Detector det_
std::pair< int, int > tileType(int layer, int ring, int phi) const
ForwardSubdetector subdet_
const HGCalDDDConstants & hdcons_
Definition: DetId.h:17
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
Detector
Definition: DetId.h:24

◆ geomDenseId2decId()

HGCalTopology::DecodedDetId HGCalTopology::geomDenseId2decId ( const uint32_t &  hi) const

Definition at line 524 of file HGCalTopology.cc.

References cellMax_, firstLay_, l1ctLayer2EG_cff::id, createfilelist::int, kHGeomHalf_, layers_, sectors_, tileTrapezoid(), totalGeomModules(), waferHexagon6(), waferMax_, and waferOff_.

524  {
526  if (hi < totalGeomModules()) {
527  id.zSide = ((int)(hi) < kHGeomHalf_ ? -1 : 1);
528  int di = ((int)(hi) % kHGeomHalf_);
529  if (waferHexagon6()) {
530  id.iSec1 = (di % sectors_);
531  di = (di - id.iSec1) / sectors_;
532  id.iLay = (di % layers_) + 1;
533  id.iType = ((di - id.iLay + 1) / layers_ == 0) ? -1 : 1;
534 #ifdef EDM_ML_DEBUG
535  edm::LogVerbatim("HGCalGeom") << "Geom Hex I/P " << hi << " O/P " << id.zSide << ":" << id.iType << ":" << id.iLay
536  << ":" << id.iSec1;
537 #endif
538  } else if (tileTrapezoid()) {
539  id.iCell1 = (di % cellMax_) + 1;
540  di = (di - id.iCell1 + 1) / cellMax_;
541  id.iSec1 = (di % sectors_) + 1;
542  di = (di - id.iSec1 + 1) / sectors_;
543  id.iLay = (di % layers_) + firstLay_;
544  id.iType = (di - id.iLay + firstLay_) / layers_;
545 #ifdef EDM_ML_DEBUG
546  edm::LogVerbatim("HGCalGeom") << "Geom Trap I/P " << hi << " O/P " << id.zSide << ":" << id.iType << ":"
547  << id.iLay << ":" << id.iSec1 << ":" << id.iCell1;
548 #endif
549  } else {
550  id.iSec2 = (di % waferMax_) - waferOff_;
551  di = (di - id.iSec2 - waferOff_) / waferMax_;
552  id.iSec1 = (di % waferMax_) - waferOff_;
553  di = (di - id.iSec1 - waferOff_) / waferMax_;
554  id.iLay = (di % layers_) + 1;
555  id.iType = (di - id.iLay + 1) / layers_;
556 #ifdef EDM_ML_DEBUG
557  edm::LogVerbatim("HGCalGeom") << "Geom Hex8 I/P " << hi << " O/P " << id.zSide << ":" << id.iType << ":"
558  << id.iLay << ":" << id.iSec1 << ":" << id.iSec2;
559 #endif
560  }
561  }
562  return id;
563 }
Log< level::Info, true > LogVerbatim
unsigned int totalGeomModules() const
Definition: HGCalTopology.h:93
bool tileTrapezoid() const
bool waferHexagon6() const
Definition: EPCuts.h:4

◆ geomMode()

HGCalGeometryMode::GeometryMode HGCalTopology::geomMode ( ) const
inline

Geometry mode.

Definition at line 79 of file HGCalTopology.h.

References mode_.

Referenced by HGCalGeometryLoader::build().

79 { return mode_; }
HGCalGeometryMode::GeometryMode mode_

◆ goEast()

DetId HGCalTopology::goEast ( const DetId id) const
inlineoverridevirtual

move the Topology east (positive ix)

Reimplemented from CaloSubdetectorTopology.

Definition at line 41 of file HGCalTopology.h.

References changeXY().

Referenced by east().

41 { return changeXY(id, +1, 0); }
DetId changeXY(const DetId &id, int nrStepsX, int nrStepsY) const
move the nagivator along x, y

◆ goNorth()

DetId HGCalTopology::goNorth ( const DetId id) const
inlineoverridevirtual

move the Topology north (increment iy)

Reimplemented from CaloSubdetectorTopology.

Definition at line 21 of file HGCalTopology.h.

References changeXY().

Referenced by north().

21 { return changeXY(id, 0, +1); }
DetId changeXY(const DetId &id, int nrStepsX, int nrStepsY) const
move the nagivator along x, y

◆ goSouth()

DetId HGCalTopology::goSouth ( const DetId id) const
inlineoverridevirtual

move the Topology south (decrement iy)

Reimplemented from CaloSubdetectorTopology.

Definition at line 31 of file HGCalTopology.h.

References changeXY().

Referenced by south().

31 { return changeXY(id, 0, -1); }
DetId changeXY(const DetId &id, int nrStepsX, int nrStepsY) const
move the nagivator along x, y

◆ goWest()

DetId HGCalTopology::goWest ( const DetId id) const
inlineoverridevirtual

move the Topology west (negative ix)

Reimplemented from CaloSubdetectorTopology.

Definition at line 51 of file HGCalTopology.h.

References changeXY().

Referenced by west().

51 { return changeXY(id, -1, 0); }
DetId changeXY(const DetId &id, int nrStepsX, int nrStepsY) const
move the nagivator along x, y

◆ isHFNose()

bool HGCalTopology::isHFNose ( ) const
inline

Definition at line 120 of file HGCalTopology.h.

References det_, DetId::Forward, HFNose, and subdet_.

Referenced by HGCalGeometryLoader::build(), HGCalGeometry::getGeometryDetId(), HGCalGeometry::getSummary(), and HGCalGeometry::newCell().

120  {
121  return (((det_ == DetId::Forward) && (subdet_ == ForwardSubdetector::HFNose)) ? true : false);
122  }
DetId::Detector det_
ForwardSubdetector subdet_

◆ maskCell()

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

Definition at line 96 of file HGCalTopology.h.

References dddConstants(), and HGCalDDDConstants::maskCell().

96 { return dddConstants().maskCell(id, corners); }
bool maskCell(const DetId &id, int corners) const
const HGCalDDDConstants & dddConstants() const
Definition: HGCalTopology.h:98

◆ neighbors()

std::vector< DetId > HGCalTopology::neighbors ( const DetId id) const

Definition at line 58 of file HGCalTopology.cc.

References addHGCSCintillatorId(), addHGCSiliconId(), HGCalCell::bottomCorner, HGCalCell::bottomLeftCorner, HGCalCell::bottomLeftEdge, HGCalCell::bottomRightCorner, HGCalCell::bottomRightEdge, HGCalDDDConstants::cellType(), HGCalCell::centralCell, decode(), HGCalDDDConstants::getTypeHex(), HGCalDDDConstants::getUVMax(), hdcons_, HGCalCell::leftEdge, HGCalDDDConstants::modifyUV(), N, HGCalCell::rightEdge, RandomServiceHelper::t1, RandomServiceHelper::t2, tileTrapezoid(), HGCalCell::topCorner, HGCalCell::topLeftCorner, HGCalCell::topLeftEdge, HGCalCell::topRightCorner, HGCalCell::topRightEdge, testProducerWithPsetDescEmpty_cfi::u1, MetAnalyzer::u2, and waferHexagon8().

58  {
59  std::vector<DetId> ids;
61  if (waferHexagon8()) {
62  int celltype = hdcons_.cellType(
63  id.iType, id.iCell1, id.iCell2, id.zSide, 1, -1); // Temporary fix - later for v17 define fwd back and orient
64 #ifdef EDM_ML_DEBUG
65  edm::LogVerbatim("HGCalGeom") << "Type:WaferU:WaferV " << id.iType << ":" << id.iCell1 << ":" << id.iCell2
66  << " CellType " << celltype;
67 #endif
68  switch (celltype) {
69  case (HGCalCell::centralCell): {
70  // cell within the wafer
71 #ifdef EDM_ML_DEBUG
72  edm::LogVerbatim("HGCalGeom") << "Cell Type 0";
73 #endif
74  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1 + 1, id.iCell2);
75  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1, id.iCell2 - 1);
76  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1 - 1, id.iCell2 - 1);
77  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1 - 1, id.iCell2);
78  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1, id.iCell2 + 1);
79  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1 + 1, id.iCell2 + 1);
80  break;
81  }
83  // bottom left edge
84  int wu1(id.iSec1), wv1(id.iSec2 - 1);
85  int t1 = hdcons_.getTypeHex(id.iLay, wu1, wv1);
86  int N1 = hdcons_.getUVMax(t1);
87  int v1 = hdcons_.modifyUV(id.iCell2, id.iType, t1);
88 #ifdef EDM_ML_DEBUG
89  edm::LogVerbatim("HGCalGeom") << "Cell Type 1 "
90  << ":" << wu1 << ":" << wv1 << ":" << t1 << ":" << N1 << ":" << v1;
91 #endif
92  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1 + 1, id.iCell2);
93  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1, id.iCell2 - 1);
94  addHGCSiliconId(ids, id.det, id.zSide, t1, id.iLay, wu1, wv1, 2 * N1 - 1, v1 + N1 - 1);
95  addHGCSiliconId(ids, id.det, id.zSide, t1, id.iLay, wu1, wv1, 2 * N1 - 1, v1 + N1);
96  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1, id.iCell2 + 1);
97  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1 + 1, id.iCell2 + 1);
98  break;
99  }
100  case (HGCalCell::leftEdge): {
101  // left edege
102  int wu1(id.iSec1 + 1), wv1(id.iSec2);
103  int t1 = hdcons_.getTypeHex(id.iLay, wu1, wv1);
104  int N1 = hdcons_.getUVMax(t1);
105  int u1 = hdcons_.modifyUV(id.iCell1, id.iType, t1);
106 #ifdef EDM_ML_DEBUG
107  edm::LogVerbatim("HGCalGeom") << "Cell Type 2 "
108  << ":" << wu1 << ":" << wv1 << ":" << t1 << ":" << N1 << ":" << u1;
109 #endif
110  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1 + 1, id.iCell2);
111  addHGCSiliconId(ids, id.det, id.zSide, t1, id.iLay, wu1, wv1, u1 + N1, 2 * N1 - 1);
112  addHGCSiliconId(ids, id.det, id.zSide, t1, id.iLay, wu1, wv1, u1 + N1 - 1, 2 * N1 - 1);
113  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1 - 1, id.iCell2);
114  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1, id.iCell2 + 1);
115  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1 + 1, id.iCell2 + 1);
116  break;
117  }
118  case (HGCalCell::topLeftEdge): {
119  // top left edge
120  int wu1(id.iSec1 + 1), wv1(id.iSec2 + 1);
121  int t1 = hdcons_.getTypeHex(id.iLay, wu1, wv1);
122  int N1 = hdcons_.getUVMax(t1);
123  int v1 = hdcons_.modifyUV(id.iCell2, id.iType, t1);
124 #ifdef EDM_ML_DEBUG
125  edm::LogVerbatim("HGCalGeom") << "Cell Type 3 "
126  << ":" << wu1 << ":" << wv1 << ":" << t1 << ":" << N1 << ":" << v1;
127 #endif
128  addHGCSiliconId(ids, id.det, id.zSide, t1, id.iLay, wu1, wv1, v1 + 1, v1 + N1);
129  addHGCSiliconId(ids, id.det, id.zSide, t1, id.iLay, wu1, wv1, v1, v1 + N1 - 1);
130  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1 - 1, id.iCell2 - 1);
131  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1 - 1, id.iCell2);
132  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1, id.iCell2 + 1);
133  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1 + 1, id.iCell2 + 1);
134  break;
135  }
136  case (HGCalCell::topRightEdge): {
137  // top right edge
138  int wu1(id.iSec1), wv1(id.iSec2 + 1);
139  int t1 = hdcons_.getTypeHex(id.iLay, wu1, wv1);
140  int N1 = hdcons_.getUVMax(t1);
141  int v1 = hdcons_.modifyUV(id.iCell2, id.iType, t1);
142 #ifdef EDM_ML_DEBUG
143  edm::LogVerbatim("HGCalGeom") << "Cell Type 4 "
144  << ":" << wu1 << ":" << wv1 << ":" << t1 << ":" << N1 << ":" << v1;
145 #endif
146  addHGCSiliconId(ids, id.det, id.zSide, t1, id.iLay, wu1, wv1, 0, v1 - N1);
147  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1, id.iCell2 - 1);
148  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1 - 1, id.iCell2 - 1);
149  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1 - 1, id.iCell2);
150  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1, id.iCell2 + 1);
151  addHGCSiliconId(ids, id.det, id.zSide, t1, id.iLay, wu1, wv1, 0, v1 - N1 + 1);
152  break;
153  }
154  case (HGCalCell::rightEdge): {
155  // right edge
156  int wu1(id.iSec1 - 1), wv1(id.iSec2);
157  int t1 = hdcons_.getTypeHex(id.iLay, wu1, wv1);
158  int N1 = hdcons_.getUVMax(t1);
159  int u1 = hdcons_.modifyUV(id.iCell1, id.iType, t1);
160 #ifdef EDM_ML_DEBUG
161  edm::LogVerbatim("HGCalGeom") << "Cell Type 5 "
162  << ":" << wu1 << ":" << wv1 << ":" << t1 << ":" << N1 << ":" << u1;
163 #endif
164  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1 + 1, id.iCell2);
165  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1, id.iCell2 - 1);
166  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1 - 1, id.iCell2 - 1);
167  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1 - 1, id.iCell2);
168  addHGCSiliconId(ids, id.det, id.zSide, t1, id.iLay, wu1, wv1, u1 - N1, 0);
169  addHGCSiliconId(ids, id.det, id.zSide, t1, id.iLay, wu1, wv1, u1 - N1 + 1, 0);
170  break;
171  }
173  // bottom right edge
174  int wu1(id.iSec1 - 1), wv1(id.iSec2 - 1);
175  int t1 = hdcons_.getTypeHex(id.iLay, wu1, wv1);
176  int N1 = hdcons_.getUVMax(t1);
177  int u1 = hdcons_.modifyUV(id.iCell1, id.iType, t1);
178 #ifdef EDM_ML_DEBUG
179  edm::LogVerbatim("HGCalGeom") << "Cell Type 6 "
180  << ":" << wu1 << ":" << wv1 << ":" << t1 << ":" << N1 << ":" << u1;
181 #endif
182  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1 + 1, id.iCell2);
183  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1, id.iCell2 - 1);
184  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1 - 1, id.iCell2 - 1);
185  addHGCSiliconId(ids, id.det, id.zSide, t1, id.iLay, wu1, wv1, u1 + N1 - 1, u1 - 1);
186  addHGCSiliconId(ids, id.det, id.zSide, t1, id.iLay, wu1, wv1, u1 + N1, u1);
187  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1 + 1, id.iCell2 + 1);
188  break;
189  }
190  case (HGCalCell::bottomCorner): {
191  // bottom corner
192  int wu1(id.iSec1), wv1(id.iSec2 - 1);
193  int t1 = hdcons_.getTypeHex(id.iLay, wu1, wv1);
194  int N1 = hdcons_.getUVMax(t1);
195  int v1 = hdcons_.modifyUV(id.iCell2, id.iType, t1);
196  int wu2(id.iSec1 - 1), wv2(id.iSec2 - 1);
197  int t2 = hdcons_.getTypeHex(id.iLay, wu2, wv2);
198  int N2 = hdcons_.getUVMax(t2);
199  int u2 = hdcons_.modifyUV(id.iCell1, id.iType, t2);
200 #ifdef EDM_ML_DEBUG
201  edm::LogVerbatim("HGCalGeom") << "Cell Type 11 "
202  << ":" << wu1 << ":" << wv1 << ":" << t1 << ":" << N1 << ":" << v1 << ":" << t2
203  << ":" << N2 << ":" << u2;
204 #endif
205  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1 + 1, id.iCell2);
206  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1, id.iCell2 - 1);
207  addHGCSiliconId(ids, id.det, id.zSide, t1, id.iLay, wu1, wv1, 2 * N1 - 1, v1 + N1 - 1);
208  addHGCSiliconId(ids, id.det, id.zSide, t1, id.iLay, wu1, wv1, 2 * N1 - 1, v1 + N1);
209  addHGCSiliconId(ids, id.det, id.zSide, t2, id.iLay, wu2, wv2, u2 + N2, u2);
210  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1 + 1, id.iCell2 + 1);
211  break;
212  }
214  // bottom left corner
215  int wu1(id.iSec1 + 1), wv1(id.iSec2);
216  int t1 = hdcons_.getTypeHex(id.iLay, wu1, wv1);
217  int N1 = hdcons_.getUVMax(t1);
218  int u1 = hdcons_.modifyUV(id.iCell1, id.iType, t1);
219  int wu2(id.iSec1), wv2(id.iSec2 - 1);
220  int t2 = hdcons_.getTypeHex(id.iLay, wu2, wv2);
221  int N2 = hdcons_.getUVMax(t2);
222  int v2 = hdcons_.modifyUV(id.iCell2, id.iType, t2);
223 #ifdef EDM_ML_DEBUG
224  edm::LogVerbatim("HGCalGeom") << "Cell Type 12 "
225  << ":" << wu1 << ":" << wv1 << ":" << t1 << ":" << N1 << ":" << u1 << ":" << t2
226  << ":" << N2 << ":" << v2;
227 #endif
228  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1 + 1, id.iCell2);
229  addHGCSiliconId(ids, id.det, id.zSide, t1, id.iLay, wu1, wv1, u1 + N1, 2 * N1 - 1);
230  addHGCSiliconId(ids, id.det, id.zSide, t2, id.iLay, wu2, wv2, 2 * N2 - 1, v2 + N2 - 1);
231  addHGCSiliconId(ids, id.det, id.zSide, t2, id.iLay, wu2, wv2, 2 * N2 - 1, v2 + N2);
232  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1, id.iCell2 + 1);
233  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1 + 1, id.iCell2 + 1);
234  break;
235  }
236  case (HGCalCell::topLeftCorner): {
237  // top left corner
238  int wu1(id.iSec1 + 1), wv1(id.iSec2 + 1);
239  int t1 = hdcons_.getTypeHex(id.iLay, wu1, wv1);
240  int N1 = hdcons_.getUVMax(t1);
241  int v1 = hdcons_.modifyUV(id.iCell2, id.iType, t1);
242  int wu2(id.iSec1 + 1), wv2(id.iSec2);
243  int t2 = hdcons_.getTypeHex(id.iLay, wu2, wv2);
244  int N2 = hdcons_.getUVMax(t2);
245  int u2 = hdcons_.modifyUV(id.iCell1, id.iType, t2);
246 #ifdef EDM_ML_DEBUG
247  edm::LogVerbatim("HGCalGeom") << "Cell Type 13 "
248  << ":" << wu1 << ":" << wv1 << ":" << t1 << ":" << N1 << ":" << v1 << ":" << t2
249  << ":" << N2 << ":" << u2;
250 #endif
251  addHGCSiliconId(ids, id.det, id.zSide, t1, id.iLay, wu1, wv1, v1 + 1, N1 + v1);
252  addHGCSiliconId(ids, id.det, id.zSide, t1, id.iLay, wu1, wv1, v1, N1 + v1 - 1);
253  addHGCSiliconId(ids, id.det, id.zSide, t2, id.iLay, wu2, wv2, u2 + N2 - 1, 2 * N2 - 1);
254  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1 - 1, id.iCell2);
255  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1, id.iCell2 + 1);
256  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1 + 1, id.iCell2 + 1);
257  break;
258  }
259  case (HGCalCell::topCorner): {
260  // top corner
261  int wu1(id.iSec1 + 1), wv1(id.iSec2 + 1);
262  int t1 = hdcons_.getTypeHex(id.iLay, wu1, wv1);
263  int N1 = hdcons_.getUVMax(t1);
264  int v1 = hdcons_.modifyUV(id.iCell2, id.iType, t1);
265  int wu2(id.iSec1), wv2(id.iSec2 + 1);
266  int t2 = hdcons_.getTypeHex(id.iLay, wu2, wv2);
267  int N2 = hdcons_.getUVMax(t2);
268  int v2 = hdcons_.modifyUV(id.iCell2, id.iType, t2);
269 #ifdef EDM_ML_DEBUG
270  edm::LogVerbatim("HGCalGeom") << "Cell Type 14 "
271  << ":" << wu1 << ":" << wv1 << ":" << t1 << ":" << N1 << ":" << v1 << ":" << t2
272  << ":" << N2 << ":" << v2;
273 #endif
274  addHGCSiliconId(ids, id.det, id.zSide, t1, id.iLay, wu1, wv1, v1 + 1, v1 + N1);
275  addHGCSiliconId(ids, id.det, id.zSide, t1, id.iLay, wu1, wv1, v1, v1 + N1 - 1);
276  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1 - 1, id.iCell2 - 1);
277  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1 - 1, id.iCell2);
278  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1, id.iCell2 + 1);
279  addHGCSiliconId(ids, id.det, id.zSide, t2, id.iLay, wu2, wv2, 0, v2 - N2 + 1);
280  break;
281  }
282  case (HGCalCell::topRightCorner): {
283  // top right corner
284  int wu1(id.iSec1), wv1(id.iSec2 + 1);
285  int t1 = hdcons_.getTypeHex(id.iLay, wu1, wv1);
286  int N1 = hdcons_.getUVMax(t1);
287  int v1 = hdcons_.modifyUV(id.iCell2, id.iType, t1);
288  int wu2(id.iSec1 - 1), wv2(id.iSec2);
289  int t2 = hdcons_.getTypeHex(id.iLay, wu2, wv2);
290  int N2 = hdcons_.getUVMax(t2);
291  int u2 = hdcons_.modifyUV(id.iCell1, id.iType, t2);
292 #ifdef EDM_ML_DEBUG
293  edm::LogVerbatim("HGCalGeom") << "Cell Type 15 "
294  << ":" << wu1 << ":" << wv1 << ":" << t1 << ":" << N1 << ":" << v1 << ":" << t2
295  << ":" << N2 << ":" << u2;
296 #endif
297  addHGCSiliconId(ids, id.det, id.zSide, t1, id.iLay, wu1, wv1, 0, v1 - N1);
298  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1, id.iCell2 - 1);
299  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1 - 1, id.iCell2 - 1);
300  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1 - 1, id.iCell2);
301  addHGCSiliconId(ids, id.det, id.zSide, t2, id.iLay, wu2, wv2, u2 - N2, 0);
302  addHGCSiliconId(ids, id.det, id.zSide, t2, id.iLay, wu2, wv2, u2 - N2 + 1, 0);
303  break;
304  }
306  // bottom right corner
307  int wu1(id.iSec1 - 1), wv1(id.iSec2 - 1);
308  int t1 = hdcons_.getTypeHex(id.iLay, wu1, wv1);
309  int N1 = hdcons_.getUVMax(t1);
310  int u1 = hdcons_.modifyUV(id.iCell1, id.iType, t1);
311  int wu2(id.iSec1 - 1), wv2(id.iSec2);
312  int t2 = hdcons_.getTypeHex(id.iLay, wu2, wv2);
313  int N2 = hdcons_.getUVMax(t2);
314  int u2 = hdcons_.modifyUV(id.iCell1, id.iType, t2);
315 #ifdef EDM_ML_DEBUG
316  edm::LogVerbatim("HGCalGeom") << "Cell Type 16 "
317  << ":" << wu1 << ":" << wv1 << ":" << t1 << ":" << N1 << ":" << u1 << ":" << t2
318  << ":" << N2 << ":" << u2;
319 #endif
320  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1 + 1, id.iCell2);
321  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1, id.iCell2 - 1);
322  addHGCSiliconId(ids, id.det, id.zSide, id.iType, id.iLay, id.iSec1, id.iSec2, id.iCell1 - 1, id.iCell2 - 1);
323  addHGCSiliconId(ids, id.det, id.zSide, t1, id.iLay, wu1, wv1, u1 + N1 - 1, u1 - 1);
324  addHGCSiliconId(ids, id.det, id.zSide, t2, id.iLay, wu2, wv2, u2 - N2, 0);
325  addHGCSiliconId(ids, id.det, id.zSide, t2, id.iLay, wu2, wv2, u2 - N2 + 1, 0);
326  break;
327  }
328  default:
329  // Not valid u, v
330  int N = hdcons_.getUVMax(id.iType);
331  edm::LogWarning("HGCalGeom") << "u:v " << id.iCell1 << ":" << id.iCell2 << " Tests " << (id.iCell1 > 2 * N - 1)
332  << ":" << (id.iCell2 > 2 * N - 1) << ":" << (id.iCell2 >= (id.iCell1 + N)) << ":"
333  << (id.iCell1 > (id.iCell2 + N)) << " ERROR";
334  }
335  } else if (tileTrapezoid()) {
336  int iphi1 = (id.iCell1 > 1) ? id.iCell1 - 1 : hdcons_.getUVMax(id.iType);
337  int iphi2 = (id.iCell1 < hdcons_.getUVMax(id.iType)) ? id.iCell1 + 1 : 1;
338  addHGCSCintillatorId(ids, id.zSide, id.iType, id.iLay, id.iSec1 - 1, id.iCell1);
339  addHGCSCintillatorId(ids, id.zSide, id.iType, id.iLay, id.iSec1 - 1, iphi1);
340  addHGCSCintillatorId(ids, id.zSide, id.iType, id.iLay, id.iSec1, iphi1);
341  addHGCSCintillatorId(ids, id.zSide, id.iType, id.iLay, id.iSec1 + 1, iphi1);
342  addHGCSCintillatorId(ids, id.zSide, id.iType, id.iLay, id.iSec1 + 1, id.iCell1);
343  addHGCSCintillatorId(ids, id.zSide, id.iType, id.iLay, id.iSec1 + 1, iphi2);
344  addHGCSCintillatorId(ids, id.zSide, id.iType, id.iLay, id.iSec1, iphi2);
345  addHGCSCintillatorId(ids, id.zSide, id.iType, id.iLay, id.iSec1 - 1, iphi2);
346  }
347  return ids;
348 }
Log< level::Info, true > LogVerbatim
static constexpr int32_t topCorner
Definition: HGCalCell.h:44
bool tileTrapezoid() const
static constexpr int32_t leftEdge
Definition: HGCalCell.h:36
int32_t cellType(int type, int waferU, int waferV, int iz, int fwdBack, int orient) const
static constexpr int32_t topLeftCorner
Definition: HGCalCell.h:43
int getTypeHex(int layer, int waferU, int waferV) const
static constexpr int32_t bottomLeftCorner
Definition: HGCalCell.h:42
int getUVMax(int type) const
DecodedDetId decode(const DetId &id) const
static constexpr int32_t topRightCorner
Definition: HGCalCell.h:45
static constexpr int32_t bottomRightEdge
Definition: HGCalCell.h:40
const HGCalDDDConstants & hdcons_
static constexpr int32_t bottomLeftEdge
Definition: HGCalCell.h:35
void addHGCSCintillatorId(std::vector< DetId > &ids, int zside, int type, int lay, int iradius, int iphi) const
add DetId of Scintillator and Silicon type if valid
static constexpr int32_t bottomRightCorner
Definition: HGCalCell.h:46
#define N
Definition: blowfish.cc:9
static constexpr int32_t centralCell
Definition: HGCalCell.h:34
static constexpr int32_t topRightEdge
Definition: HGCalCell.h:38
bool waferHexagon8() const
void addHGCSiliconId(std::vector< DetId > &ids, int det, int zside, int type, int lay, int waferU, int waferV, int cellU, int cellV) const
static constexpr int32_t rightEdge
Definition: HGCalCell.h:39
static constexpr int32_t bottomCorner
Definition: HGCalCell.h:41
int modifyUV(int uv, int type1, int type2) const
Log< level::Warning, false > LogWarning
static constexpr int32_t topLeftEdge
Definition: HGCalCell.h:37

◆ north()

std::vector<DetId> HGCalTopology::north ( const DetId id) const
inlineoverridevirtual

Get the neighbors of the given cell in north direction

Implements CaloSubdetectorTopology.

Definition at line 22 of file HGCalTopology.h.

References goNorth().

22  {
23  DetId nextId = goNorth(id);
24  std::vector<DetId> vNeighborsDetId;
25  if (!(nextId == DetId(0)))
26  vNeighborsDetId.emplace_back(nextId);
27  return vNeighborsDetId;
28  }
DetId goNorth(const DetId &id) const override
move the Topology north (increment iy)
Definition: HGCalTopology.h:21
Definition: DetId.h:17

◆ offsetBy()

DetId HGCalTopology::offsetBy ( const DetId  startId,
int  nrStepsX,
int  nrStepsY 
) const

returns a new DetId offset by nrStepsX and nrStepsY (can be negative), returns DetId(0) if invalid

Definition at line 505 of file HGCalTopology.cc.

References changeXY(), DetId::det(), DetId::Forward, l1ctLayer2EG_cff::id, createfilelist::int, subdet_, DetId::subdetId(), and valid().

505  {
506  if (startId.det() == DetId::Forward && startId.subdetId() == (int)(subdet_)) {
507  DetId id = changeXY(startId, nrStepsX, nrStepsY);
508  if (valid(id))
509  return id;
510  }
511  return DetId(0);
512 }
DetId changeXY(const DetId &id, int nrStepsX, int nrStepsY) const
move the nagivator along x, y
bool valid(const DetId &id) const override
Is this a valid cell id.
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:46
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
ForwardSubdetector subdet_
Definition: DetId.h:17

◆ south()

std::vector<DetId> HGCalTopology::south ( const DetId id) const
inlineoverridevirtual

Get the neighbors of the given cell in south direction

Implements CaloSubdetectorTopology.

Definition at line 32 of file HGCalTopology.h.

References goSouth().

32  {
33  DetId nextId = goSouth(id);
34  std::vector<DetId> vNeighborsDetId;
35  if (!(nextId == DetId(0)))
36  vNeighborsDetId.emplace_back(nextId);
37  return vNeighborsDetId;
38  }
Definition: DetId.h:17
DetId goSouth(const DetId &id) const override
move the Topology south (decrement iy)
Definition: HGCalTopology.h:31

◆ subDetector()

ForwardSubdetector HGCalTopology::subDetector ( ) const
inline

Definition at line 118 of file HGCalTopology.h.

References subdet_.

Referenced by HGCalGeometryLoader::build().

118 { return subdet_; }
ForwardSubdetector subdet_

◆ switchZSide()

DetId HGCalTopology::switchZSide ( const DetId  startId) const

Definition at line 514 of file HGCalTopology.cc.

References decode(), encode(), l1ctLayer2EG_cff::id, valid(), and HGCalTopology::DecodedDetId::zSide.

514  {
515  HGCalTopology::DecodedDetId id_ = decode(startId);
516  id_.zSide = -id_.zSide;
517  DetId id = encode(id_);
518  if (valid(id))
519  return id;
520  else
521  return DetId(0);
522 }
bool valid(const DetId &id) const override
Is this a valid cell id.
DecodedDetId decode(const DetId &id) const
Definition: DetId.h:17
DetId encode(const DecodedDetId &id_) const

◆ tileTrapezoid()

bool HGCalTopology::tileTrapezoid ( ) const
inline

◆ totalGeomModules()

unsigned int HGCalTopology::totalGeomModules ( ) const
inline

◆ totalModules()

unsigned int HGCalTopology::totalModules ( ) const
inline

Definition at line 92 of file HGCalTopology.h.

References kSizeForDenseIndexing.

Referenced by HGCalGeometry::HGCalGeometry().

92 { return kSizeForDenseIndexing; }
unsigned int kSizeForDenseIndexing

◆ up()

std::vector<DetId> HGCalTopology::up ( const DetId id) const
inlineoverridevirtual

Get the neighbors of the given cell in up direction (outward)

Implements CaloSubdetectorTopology.

Definition at line 60 of file HGCalTopology.h.

References changeZ().

60  {
61  DetId nextId = changeZ(id, +1);
62  std::vector<DetId> vNeighborsDetId;
63  if (!(nextId == DetId(0)))
64  vNeighborsDetId.emplace_back(nextId);
65  return vNeighborsDetId;
66  }
Definition: DetId.h:17
DetId changeZ(const DetId &id, int nrStepsZ) const
move the nagivator along z

◆ valid() [1/2]

bool HGCalTopology::valid ( const DetId id) const
overridevirtual

Is this a valid cell id.

Reimplemented from CaloSubdetectorTopology.

Definition at line 465 of file HGCalTopology.cc.

References cells_, decode(), DetId::det(), det_, RemoveAddSevLevel::flag, hdcons_, createfilelist::int, HGCalDDDConstants::isValidHex(), HGCalDDDConstants::isValidHex8(), HGCalDDDConstants::isValidTrap(), layers_, sectors_, subdet_, DetId::subdetId(), tileTrapezoid(), and waferHexagon6().

Referenced by HGCDigitizer::checkPosition(), HGCalGeometry::newCell(), offsetBy(), switchZSide(), valid(), HGCalTriggerGeometryV9Imp2::validCell(), HGCalTriggerGeometryV9Imp3::validCell(), HGCalTriggerGeometryV9Imp2::validCellId(), HGCalTriggerGeometryV9Imp3::validCellId(), and validModule().

465  {
467  bool flag;
468  if (waferHexagon6()) {
469  flag = (idin.det() == det_ && idin.subdetId() == (int)(subdet_) && id.iCell1 >= 0 && id.iCell1 < cells_ &&
470  id.iLay > 0 && id.iLay <= layers_ && id.iSec1 >= 0 && id.iSec1 <= sectors_);
471  if (flag)
472  flag = hdcons_.isValidHex(id.iLay, id.iSec1, id.iCell1, true);
473  } else if (tileTrapezoid()) {
474  flag = ((idin.det() == det_) && hdcons_.isValidTrap(id.zSide, id.iLay, id.iSec1, id.iCell1));
475  } else {
476  flag = ((idin.det() == det_) && hdcons_.isValidHex8(id.iLay, id.iSec1, id.iSec2, id.iCell1, id.iCell2));
477  }
478  return flag;
479 }
bool tileTrapezoid() const
bool isValidHex(int lay, int mod, int cell, bool reco) const
bool isValidHex8(int lay, int waferU, int waferV, bool fullAndPart=false) const
bool waferHexagon6() const
DetId::Detector det_
DecodedDetId decode(const DetId &id) const
bool isValidTrap(int zside, int lay, int ieta, int iphi) const
ForwardSubdetector subdet_
const HGCalDDDConstants & hdcons_

◆ valid() [2/2]

bool HGCalTopology::valid ( const DetId id,
int  cornerMin 
) const

Definition at line 481 of file HGCalTopology.cc.

References decode(), DetId::det(), det_, RemoveAddSevLevel::flag, hdcons_, HGCalDDDConstants::isValidHex8(), gpuClustering::pixelStatus::mask, HGCalDDDConstants::maskCell(), valid(), HGCalTypes::WaferCornerMin, and waferHexagon8().

481  {
482  if (waferHexagon8()) {
484  bool mask = (cornerMin < HGCalTypes::WaferCornerMin) ? false : hdcons_.maskCell(idin, cornerMin);
485  bool flag = ((idin.det() == det_) &&
487  id.iLay, id.iSec1, id.iSec2, id.iCell1, id.iCell2, (cornerMin >= HGCalTypes::WaferCornerMin)));
488  return (flag && (!mask));
489  } else {
490  return valid(idin);
491  }
492 }
bool maskCell(const DetId &id, int corners) const
bool valid(const DetId &id) const override
Is this a valid cell id.
bool isValidHex8(int lay, int waferU, int waferV, bool fullAndPart=false) const
DetId::Detector det_
constexpr uint32_t mask
Definition: gpuClustering.h:24
DecodedDetId decode(const DetId &id) const
const HGCalDDDConstants & hdcons_
bool waferHexagon8() const
static constexpr int32_t WaferCornerMin
Definition: HGCalTypes.h:74

◆ validHashIndex()

bool HGCalTopology::validHashIndex ( uint32_t  ix) const
inline

Definition at line 89 of file HGCalTopology.h.

References kSizeForDenseIndexing.

Referenced by denseId2detId().

89 { return (ix < kSizeForDenseIndexing); }
unsigned int kSizeForDenseIndexing

◆ validModule()

bool HGCalTopology::validModule ( const DetId id,
int  cornerMin 
) const

Definition at line 494 of file HGCalTopology.cc.

References decode(), DetId::det(), det_, hdcons_, DetId::HGCalEE, DetId::HGCalHSi, HGCalDDDConstants::isValidHex8(), valid(), and HGCalTypes::WaferCornerMin.

494  {
495  if (idin.det() != det_) {
496  return false;
497  } else if ((idin.det() == DetId::HGCalEE) || (idin.det() == DetId::HGCalHSi)) {
499  return hdcons_.isValidHex8(id.iLay, id.iSec1, id.iSec2, (cornerMin >= HGCalTypes::WaferCornerMin));
500  } else {
501  return valid(idin);
502  }
503 }
bool valid(const DetId &id) const override
Is this a valid cell id.
bool isValidHex8(int lay, int waferU, int waferV, bool fullAndPart=false) const
DetId::Detector det_
DecodedDetId decode(const DetId &id) const
const HGCalDDDConstants & hdcons_
static constexpr int32_t WaferCornerMin
Definition: HGCalTypes.h:74

◆ waferHexagon6()

bool HGCalTopology::waferHexagon6 ( ) const
inline

◆ waferHexagon8()

bool HGCalTopology::waferHexagon8 ( ) const
inline

Definition at line 126 of file HGCalTopology.h.

References hdcons_, and HGCalDDDConstants::waferHexagon8().

Referenced by HGCalGeometry::getClosestCellHex(), neighbors(), and valid().

126 { return hdcons_.waferHexagon8(); }
bool waferHexagon8() const
const HGCalDDDConstants & hdcons_

◆ west()

std::vector<DetId> HGCalTopology::west ( const DetId id) const
inlineoverridevirtual

Get the neighbors of the given cell in west direction

Implements CaloSubdetectorTopology.

Definition at line 52 of file HGCalTopology.h.

References goWest().

52  {
53  DetId nextId = goWest(id);
54  std::vector<DetId> vNeighborsDetId;
55  if (!(nextId == DetId(0)))
56  vNeighborsDetId.emplace_back(nextId);
57  return vNeighborsDetId;
58  }
Definition: DetId.h:17
DetId goWest(const DetId &id) const override
move the Topology west (negative ix)
Definition: HGCalTopology.h:51

Member Data Documentation

◆ cellMax_

int HGCalTopology::cellMax_
private

◆ cells_

int HGCalTopology::cells_
private

Definition at line 146 of file HGCalTopology.h.

Referenced by HGCalTopology(), and valid().

◆ det_

DetId::Detector HGCalTopology::det_
private

Definition at line 144 of file HGCalTopology.h.

Referenced by decode(), detector(), encode(), HGCalTopology(), isHFNose(), valid(), and validModule().

◆ firstLay_

int HGCalTopology::firstLay_
private

◆ hdcons_

const HGCalDDDConstants& HGCalTopology::hdcons_
private

◆ kHGeomHalf_

int HGCalTopology::kHGeomHalf_
private

◆ kHGhalf_

int HGCalTopology::kHGhalf_
private

Definition at line 148 of file HGCalTopology.h.

Referenced by HGCalTopology().

◆ kHGhalfType_

int HGCalTopology::kHGhalfType_
private

Definition at line 148 of file HGCalTopology.h.

Referenced by denseId2detId(), detId2denseId(), and HGCalTopology().

◆ kSizeForDenseIndexing

unsigned int HGCalTopology::kSizeForDenseIndexing
private

Definition at line 149 of file HGCalTopology.h.

Referenced by HGCalTopology(), totalModules(), and validHashIndex().

◆ layers_

int HGCalTopology::layers_
private

◆ mode_

HGCalGeometryMode::GeometryMode HGCalTopology::mode_
private

Definition at line 142 of file HGCalTopology.h.

Referenced by geomMode(), and HGCalTopology().

◆ sectors_

int HGCalTopology::sectors_
private

◆ subdet_

ForwardSubdetector HGCalTopology::subdet_
private

Definition at line 145 of file HGCalTopology.h.

Referenced by decode(), encode(), HGCalTopology(), isHFNose(), offsetBy(), subDetector(), and valid().

◆ subSectors_

const int HGCalTopology::subSectors_ = 2
static

Use subSector in square mode as wafer type in hexagon mode.

Definition at line 106 of file HGCalTopology.h.

◆ types_

int HGCalTopology::types_
private

Definition at line 146 of file HGCalTopology.h.

Referenced by denseId2detId(), detId2denseId(), and HGCalTopology().

◆ waferMax_

int HGCalTopology::waferMax_
private

◆ waferOff_

int HGCalTopology::waferOff_
private