CMS 3D CMS Logo

HGCSiliconDetId.h
Go to the documentation of this file.
1 #ifndef DataFormats_ForwardDetId_HGCSiliconDetId_H
2 #define DataFormats_ForwardDetId_HGCSiliconDetId_H 1
3 
4 #include <iosfwd>
8 
9 /* \brief description of the bit assigment
10  [0:4] u-coordinate of the cell (measured from the lower left
11  [5:9] v-coordinate of the cell corner of the wafer)
12  [10:13] abs(u) of the wafer (u-axis points along -x axis)
13  [14:14] sign of u (0:+u; 1:-u) (u=0 is at the center of beam line)
14  [15:18] abs(v) of the wafer (v-axis points 60-degree wrt x-axis)
15  [19:19] sign of v (0:+v; 1:-v) (v=0 is at the center of beam line)
16  [20:24] layer number
17  [25:25] z-side (0 for +z; 1 for -z)
18  [26:27] Type (0 fine divisions of wafer with 120 mum thick silicon
19  1 coarse divisions of wafer with 200 mum thick silicon
20  2 coarse divisions of wafer with 300 mum thick silicon)
21  [28:31] Detector type (HGCalEE or HGCalHSi)
22 */
23 class HGCSiliconDetId : public DetId {
24 public:
26  static const int HGCalFineN = 12;
27  static const int HGCalCoarseN = 8;
28  static const int HGCalFineTrigger = 3;
29  static const int HGCalCoarseTrigger = 2;
30 
32  constexpr HGCSiliconDetId() : DetId() {}
34  constexpr HGCSiliconDetId(uint32_t rawid) : DetId(rawid) {}
36  constexpr HGCSiliconDetId(
37  DetId::Detector det, int zp, int type, int layer, int waferU, int waferV, int cellU, int cellV)
38  : DetId(det, ForwardEmpty) {
39  int waferUabs(std::abs(waferU)), waferVabs(std::abs(waferV));
40  int waferUsign = (waferU >= 0) ? 0 : 1;
41  int waferVsign = (waferV >= 0) ? 0 : 1;
42  int zside = (zp < 0) ? 1 : 0;
44  ((waferUabs & kHGCalWaferUMask) << kHGCalWaferUOffset) |
45  ((waferUsign & kHGCalWaferUSignMask) << kHGCalWaferUSignOffset) |
46  ((waferVabs & kHGCalWaferVMask) << kHGCalWaferVOffset) |
47  ((waferVsign & kHGCalWaferVSignMask) << kHGCalWaferVSignOffset) |
50  }
51 
53  constexpr HGCSiliconDetId(const DetId& gen) {
54  if (!gen.null()) {
55  if ((gen.det() != HGCalEE) && (gen.det() != HGCalHSi)) {
56  throw cms::Exception("Invalid DetId")
57  << "Cannot initialize HGCSiliconDetId from " << std::hex << gen.rawId() << std::dec;
58  }
59  }
60  id_ = gen.rawId();
61  }
62 
64  constexpr HGCSiliconDetId& operator=(const DetId& gen) {
65  if (!gen.null()) {
66  if ((gen.det() != HGCalEE) && (gen.det() != HGCalHSi)) {
67  throw cms::Exception("Invalid DetId")
68  << "Cannot assign HGCSiliconDetId from " << std::hex << gen.rawId() << std::dec;
69  }
70  }
71  id_ = gen.rawId();
72  return (*this);
73  }
74 
76  constexpr HGCSiliconDetId geometryCell() const {
77  return HGCSiliconDetId(det(), zside(), 0, layer(), waferU(), waferV(), 0, 0);
78  }
79  constexpr HGCSiliconDetId moduleId() const {
80  return HGCSiliconDetId(det(), zside(), type(), layer(), waferU(), waferV(), 0, 0);
81  }
82 
84  constexpr DetId::Detector subdet() const { return det(); }
85 
87  constexpr int type() const { return (id_ >> kHGCalTypeOffset) & kHGCalTypeMask; }
88 
90  constexpr int zside() const { return (((id_ >> kHGCalZsideOffset) & kHGCalZsideMask) ? -1 : 1); }
91 
93  constexpr int layer() const { return (id_ >> kHGCalLayerOffset) & kHGCalLayerMask; }
94 
96  constexpr int cellU() const { return (id_ >> kHGCalCellUOffset) & kHGCalCellUMask; }
97  constexpr int cellV() const { return (id_ >> kHGCalCellVOffset) & kHGCalCellVMask; }
98  constexpr std::pair<int, int> cellUV() const { return std::pair<int, int>(cellU(), cellV()); }
99  constexpr int cellX() const {
100  int N = (type() == HGCalFine) ? HGCalFineN : HGCalCoarseN;
101  return (3 * (cellV() - N) + 2);
102  }
103  constexpr int cellY() const {
104  int N = (type() == HGCalFine) ? HGCalFineN : HGCalCoarseN;
105  return (2 * cellU() - (N + cellV()));
106  }
107  constexpr std::pair<int, int> cellXY() const { return std::pair<int, int>(cellX(), cellY()); }
108 
110  constexpr int waferUAbs() const { return (id_ >> kHGCalWaferUOffset) & kHGCalWaferUMask; }
111  constexpr int waferVAbs() const { return (id_ >> kHGCalWaferVOffset) & kHGCalWaferVMask; }
112  constexpr int waferU() const {
114  }
115  constexpr int waferV() const {
117  }
118  constexpr std::pair<int, int> waferUV() const { return std::pair<int, int>(waferU(), waferV()); }
119  constexpr int waferX() const { return (-2 * waferU() + waferV()); }
120  constexpr int waferY() const { return (2 * waferV()); }
121  constexpr std::pair<int, int> waferXY() const { return std::pair<int, int>(waferX(), waferY()); }
122  constexpr void unpack(int& ty, int& zs, int& ly, int& wU, int& wV, int& cU, int& cV) const {
123  ty = type();
124  zs = zside();
125  ly = layer();
126  wU = waferU();
127  wV = waferV();
128  cU = cellU();
129  cV = cellV();
130  }
131 
132  // get trigger cell u,v
133  constexpr int triggerCellU() const {
134  int N = (type() == HGCalFine) ? HGCalFineN : HGCalCoarseN;
136  return (cellU() >= N && cellV() >= N)
137  ? cellU() / NT
138  : ((cellU() < N && cellU() <= cellV()) ? cellU() / NT : (1 + (cellU() - (cellV() % NT + 1)) / NT));
139  }
140  constexpr int triggerCellV() const {
141  int N = (type() == HGCalFine) ? HGCalFineN : HGCalCoarseN;
143  return (cellU() >= N && cellV() >= N)
144  ? cellV() / NT
145  : ((cellU() < N && cellU() <= cellV()) ? ((cellV() - cellU()) / NT + cellU() / NT) : cellV() / NT);
146  }
147  constexpr std::pair<int, int> triggerCellUV() const { return std::pair<int, int>(triggerCellU(), triggerCellV()); }
148 
150  constexpr bool isEE() const { return (det() == HGCalEE); }
151  constexpr bool isHE() const { return (det() == HGCalHSi); }
152  constexpr bool isForward() const { return true; }
153 
155 
156 public:
157  static constexpr uint32_t kHGCalCellUOffset = 0;
158  static constexpr uint32_t kHGCalCellUMask = 0x1F;
159  static constexpr uint32_t kHGCalCellVOffset = 5;
160  static constexpr uint32_t kHGCalCellVMask = 0x1F;
161  static constexpr uint32_t kHGCalWaferUOffset = 10;
162  static constexpr uint32_t kHGCalWaferUMask = 0xF;
163  static constexpr uint32_t kHGCalWaferUSignOffset = 14;
164  static constexpr uint32_t kHGCalWaferUSignMask = 0x1;
165  static constexpr uint32_t kHGCalWaferVOffset = 15;
166  static constexpr uint32_t kHGCalWaferVMask = 0xF;
167  static constexpr uint32_t kHGCalWaferVSignOffset = 19;
168  static constexpr uint32_t kHGCalWaferVSignMask = 0x1;
169  static constexpr uint32_t kHGCalLayerOffset = 20;
170  static constexpr uint32_t kHGCalLayerMask = 0x1F;
171  static constexpr uint32_t kHGCalZsideOffset = 25;
172  static constexpr uint32_t kHGCalZsideMask = 0x1;
173  static constexpr uint32_t kHGCalTypeOffset = 26;
174  static constexpr uint32_t kHGCalTypeMask = 0x3;
175 };
176 
177 std::ostream& operator<<(std::ostream&, const HGCSiliconDetId& id);
178 
179 #endif
static constexpr uint32_t kHGCalCellUOffset
constexpr void unpack(int &ty, int &zs, int &ly, int &wU, int &wV, int &cU, int &cV) const
constexpr int cellV() const
constexpr int cellX() const
constexpr int waferUAbs() const
get the wafer #&#39;s in u,v or in x,y
static constexpr uint32_t kHGCalCellVMask
constexpr int waferY() const
constexpr int waferVAbs() const
static constexpr uint32_t kHGCalTypeOffset
constexpr HGCSiliconDetId(uint32_t rawid)
constexpr HGCSiliconDetId & operator=(const DetId &gen)
static const int HGCalFineN
constexpr bool isHE() const
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:46
constexpr std::pair< int, int > waferUV() const
constexpr DetId::Detector subdet() const
get the subdetector
static constexpr uint32_t kHGCalWaferUMask
constexpr int zside() const
get the z-side of the cell (1/-1)
static const int HGCalCoarseTrigger
constexpr std::pair< int, int > cellXY() const
static const HGCSiliconDetId Undefined
static constexpr uint32_t kHGCalWaferVOffset
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::ostream & operator<<(std::ostream &, const HGCSiliconDetId &id)
static constexpr uint32_t kHGCalWaferUOffset
constexpr int layer() const
get the layer #
static const int HGCalCoarseN
constexpr int cellU() const
get the cell #&#39;s in u,v or in x,y
static constexpr uint32_t kHGCalWaferVSignMask
constexpr HGCSiliconDetId geometryCell() const
constexpr int type() const
get the type
static constexpr uint32_t kHGCalTypeMask
static constexpr uint32_t kHGCalLayerOffset
Definition: DetId.h:17
static constexpr uint32_t kHGCalCellUMask
#define N
Definition: blowfish.cc:9
static constexpr uint32_t kHGCalZsideOffset
constexpr int waferX() const
constexpr HGCSiliconDetId moduleId() const
constexpr bool isForward() const
Detector
Definition: DetId.h:24
uint32_t id_
Definition: DetId.h:69
constexpr HGCSiliconDetId(DetId::Detector det, int zp, int type, int layer, int waferU, int waferV, int cellU, int cellV)
static constexpr uint32_t kHGCalLayerMask
constexpr int waferU() const
constexpr HGCSiliconDetId()
static constexpr uint32_t kHGCalWaferUSignMask
static constexpr uint32_t kHGCalWaferVMask
static constexpr uint32_t kHGCalWaferUSignOffset
static constexpr uint32_t kHGCalWaferVSignOffset
constexpr int triggerCellV() const
constexpr int triggerCellU() const
constexpr std::pair< int, int > triggerCellUV() const
constexpr int cellY() const
constexpr std::pair< int, int > cellUV() const
static constexpr uint32_t kHGCalZsideMask
static const int HGCalFineTrigger
constexpr std::pair< int, int > waferXY() const
constexpr bool isEE() const
consistency check : no bits left => no overhead
constexpr int waferV() const
static constexpr uint32_t kHGCalCellVOffset
constexpr HGCSiliconDetId(const DetId &gen)