CMS 3D CMS Logo

HGCalTBTopology.cc
Go to the documentation of this file.
5 
6 //#define EDM_ML_DEBUG
7 
8 HGCalTBTopology::HGCalTBTopology(const HGCalTBDDDConstants* hdcons, int det) : hdcons_(hdcons) {
10  layers_ = hdcons_->layers(true);
11  cells_ = hdcons_->maxCells(true);
12  mode_ = hdcons_->geomMode();
16  subdet_ = static_cast<ForwardSubdetector>(det);
18  types_ = 2;
20  kSizeForDenseIndexing = static_cast<unsigned int>(2 * kHGhalf_);
21 #ifdef EDM_ML_DEBUG
22  edm::LogVerbatim("HGCalGeom") << "HGCalTBTopology initialized for detector " << det << ":" << det_ << ":" << subdet_
23  << " having " << sectors_ << " Sectors, " << layers_ << " Layers from " << firstLay_
24  << ", " << cells_ << " cells and total channels " << kSizeForDenseIndexing << ":"
25  << (2 * kHGeomHalf_);
26 #endif
27 }
28 
30 
31 std::vector<DetId> HGCalTBTopology::neighbors(DetId idin) const {
32  std::vector<DetId> ids = north(idin);
33  for (const auto& id : south(idin))
34  ids.emplace_back(id);
35  for (const auto& id : east(idin))
36  ids.emplace_back(id);
37  for (const auto& id : west(idin))
38  ids.emplace_back(id);
39  return ids;
40 }
41 
42 unsigned int HGCalTBTopology::allGeomModules() const { return (static_cast<unsigned int>(2 * hdcons_->wafers())); }
43 
46  if (validHashIndex(hi)) {
47  id.zSide = ((int)(hi) < kHGhalfType_ ? -1 : 1);
48  int di = ((int)(hi) % kHGhalfType_);
49  int type = (di % types_);
50  id.iType = (type == 0 ? -1 : 1);
51  id.iSec1 = (((di - type) / types_) % sectors_);
52  id.iLay = (((((di - type) / types_) - id.iSec1 + 1) / sectors_) % layers_ + 1);
53  id.iCell1 = (((((di - type) / types_) - id.iSec1 + 1) / sectors_ - id.iLay + 1) / layers_ + 1);
54 #ifdef EDM_ML_DEBUG
55  edm::LogVerbatim("HGCalGeom") << "Input Hex " << hi << " o/p " << id.zSide << ":" << id.iLay << ":" << id.iType
56  << ":" << id.iSec1 << ":" << id.iCell1;
57 #endif
58  }
59  return encode(id);
60 }
61 
62 uint32_t HGCalTBTopology::detId2denseId(const DetId& idin) const {
64  int type = (id.iType > 0) ? 1 : 0;
65  uint32_t idx =
66  static_cast<uint32_t>(((id.zSide > 0) ? kHGhalfType_ : 0) +
67  ((((id.iCell1 - 1) * layers_ + id.iLay - 1) * sectors_ + id.iSec1) * types_ + type));
68 #ifdef EDM_ML_DEBUG
69  edm::LogVerbatim("HGCalGeom") << "Input Hex " << id.zSide << ":" << id.iLay << ":" << id.iSec1 << ":" << id.iCell1
70  << ":" << id.iType << " Constants " << kHGeomHalf_ << ":" << layers_ << ":" << sectors_
71  << ":" << types_ << " o/p " << idx;
72 #endif
73  return idx;
74 }
75 
76 uint32_t HGCalTBTopology::detId2denseGeomId(const DetId& idin) const {
78  uint32_t idx;
79  idx = (uint32_t)(((id.zSide > 0) ? kHGeomHalf_ : 0) + (id.iLay - 1) * sectors_ + id.iSec1);
80 #ifdef EDM_ML_DEBUG
81  edm::LogVerbatim("HGCalGeom") << "Geom Hex I/P " << id.zSide << ":" << id.iLay << ":" << id.iSec1 << ":" << id.iType
82  << " Constants " << kHGeomHalf_ << ":" << layers_ << ":" << sectors_ << " o/p " << idx;
83 #endif
84  return idx;
85 }
86 
87 bool HGCalTBTopology::valid(const DetId& idin) const {
89  bool flag;
90  flag = (idin.det() == det_ && idin.subdetId() == (int)(subdet_) && id.iCell1 >= 0 && id.iCell1 < cells_ &&
91  id.iLay > 0 && id.iLay <= layers_ && id.iSec1 >= 0 && id.iSec1 <= sectors_);
92  if (flag)
93  flag = hdcons_->isValidHex(id.iLay, id.iSec1, id.iCell1, true);
94  return flag;
95 }
96 
97 DetId HGCalTBTopology::offsetBy(const DetId startId, int nrStepsX, int nrStepsY) const {
98  if (startId.det() == DetId::Forward && startId.subdetId() == static_cast<int>(subdet_)) {
99  DetId id = changeXY(startId, nrStepsX, nrStepsY);
100  if (valid(id))
101  return id;
102  }
103  return DetId(0);
104 }
105 
107  HGCalTBTopology::DecodedDetId id_ = decode(startId);
108  id_.zSide = -id_.zSide;
109  DetId id = encode(id_);
110  if (valid(id))
111  return id;
112  else
113  return DetId(0);
114 }
115 
118  if (hi < totalGeomModules()) {
119  id.zSide = ((int)(hi) < kHGeomHalf_ ? -1 : 1);
120  int di = ((int)(hi) % kHGeomHalf_);
121  id.iSec1 = (di % sectors_);
122  di = (di - id.iSec1) / sectors_;
123  id.iLay = (di % layers_) + 1;
124  id.iType = ((di - id.iLay + 1) / layers_ == 0) ? -1 : 1;
125 #ifdef EDM_ML_DEBUG
126  edm::LogVerbatim("HGCalGeom") << "Geom Hex I/P " << hi << " O/P " << id.zSide << ":" << id.iType << ":" << id.iLay
127  << ":" << id.iSec1;
128 #endif
129  }
130  return id;
131 }
132 
135  HGCalDetId id(startId);
136  idx.iCell1 = id.cell();
137  idx.iCell2 = 0;
138  idx.iLay = id.layer();
139  idx.iSec1 = id.wafer();
140  idx.iSec2 = 0;
141  idx.iType = id.waferType();
142  idx.zSide = id.zside();
143  idx.det = id.subdetId();
144  return idx;
145 }
146 
148  DetId id;
149 #ifdef EDM_ML_DEBUG
150  edm::LogVerbatim("HGCalGeomX") << "Encode " << idx.det << ":" << idx.zSide << ":" << idx.iType << ":" << idx.iLay
151  << ":" << idx.iSec1 << ":" << idx.iSec2 << ":" << idx.iCell1 << ":" << idx.iCell2;
152 #endif
153  id = HGCalDetId((ForwardSubdetector)(idx.det), idx.zSide, idx.iLay, ((idx.iType > 0) ? 1 : 0), idx.iSec1, idx.iCell1)
154  .rawId();
155  return id;
156 }
157 
158 DetId HGCalTBTopology::changeXY(const DetId& id, int nrStepsX, int nrStepsY) const { return DetId(); }
159 
160 DetId HGCalTBTopology::changeZ(const DetId& id, int nrStepsZ) const { return DetId(); }
161 
163 
unsigned int totalGeomModules() const
Log< level::Info, true > LogVerbatim
unsigned int kSizeForDenseIndexing
DetId::Detector det_
HGCalTBTopology(const HGCalTBDDDConstants *hdcons, int subdet)
create a new Topology
DetId changeXY(const DetId &id, int nrStepsX, int nrStepsY) const
move the nagivator along x, y
DetId encode(const DecodedDetId &id_) const
std::vector< DetId > north(const DetId &id) const override
ForwardSubdetector
std::vector< DetId > neighbors(DetId id) const
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:46
DecodedDetId geomDenseId2decId(const uint32_t &hi) const
Definition: EPCuts.h:4
std::vector< DetId > west(const DetId &id) const override
uint32_t detId2denseId(const DetId &id) const override
return a linear packed id
~HGCalTBTopology() override
default destructor
DetId switchZSide(const DetId startId) const
DetId denseId2detId(uint32_t denseId) const override
Dense indexing.
bool validHashIndex(uint32_t ix) const
std::vector< DetId > east(const DetId &id) const override
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
DetId offsetBy(const DetId startId, int nrStepsX, int nrStepsY) const
const HGCalTBDDDConstants * hdcons_
Definition: DetId.h:17
DecodedDetId decode(const DetId &id) const
ForwardSubdetector subdet_
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
bool valid(const DetId &id) const override
Is this a valid cell id.
virtual uint32_t detId2denseGeomId(const DetId &id) const
std::vector< DetId > south(const DetId &id) const override
unsigned int allGeomModules() const
HGCalGeometryMode::GeometryMode mode_
unsigned int layers(bool reco) const
#define TYPELOOKUP_DATA_REG(_dataclass_)
Definition: typelookup.h:102
HGCalGeometryMode::GeometryMode geomMode() const
bool isValidHex(int lay, int mod, int cell, bool reco) const
int maxCells(bool reco) const
DetId changeZ(const DetId &id, int nrStepsZ) const
move the nagivator along z