CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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>
7 
8 /* \brief description of the bit assigment
9  [0:4] u-coordinate of the cell (measured from the lower left
10  [5:9] v-coordinate of the cell corner of the wafer)
11  [10:13] abs(u) of the wafer (u-axis points along -x axis)
12  [14:14] sign of u (0:+u; 1:-u) (u=0 is at the center of beam line)
13  [15:18] abs(v) of the wafer (v-axis points 60-degree wrt x-axis)
14  [19:19] sign of v (0:+v; 1:-v) (v=0 is at the center of beam line)
15  [20:24] layer number
16  [25:25] z-side (0 for +z; 1 for -z)
17  [26:27] Type (0 fine divisions of wafer with 120 mum thick silicon
18  1 coarse divisions of wafer with 200 mum thick silicon
19  2 coarse divisions of wafer with 300 mum thick silicon)
20  [28:31] Detector type (HGCalEE or HGCalHSi)
21 */
22 class HGCSiliconDetId : public DetId {
23 public:
25  static const int HGCalFineN = 12;
26  static const int HGCalCoarseN = 8;
27  static const int HGCalFineTrigger = 3;
28  static const int HGCalCoarseTrigger = 2;
29 
33  HGCSiliconDetId(uint32_t rawid);
35  HGCSiliconDetId(DetId::Detector det, int zp, int type, int layer, int waferU, int waferV, int cellU, int cellV);
37  HGCSiliconDetId(const DetId& id);
39  HGCSiliconDetId& operator=(const DetId& id);
40 
42  HGCSiliconDetId geometryCell() const { return HGCSiliconDetId(det(), zside(), 0, layer(), waferU(), waferV(), 0, 0); }
44  return HGCSiliconDetId(det(), zside(), type(), layer(), waferU(), waferV(), 0, 0);
45  }
46 
48  DetId::Detector subdet() const { return det(); }
49 
51  int type() const { return (id_ >> kHGCalTypeOffset) & kHGCalTypeMask; }
52 
54  int zside() const { return (((id_ >> kHGCalZsideOffset) & kHGCalZsideMask) ? -1 : 1); }
55 
57  int layer() const { return (id_ >> kHGCalLayerOffset) & kHGCalLayerMask; }
58 
60  int cellU() const { return (id_ >> kHGCalCellUOffset) & kHGCalCellUMask; }
61  int cellV() const { return (id_ >> kHGCalCellVOffset) & kHGCalCellVMask; }
62  std::pair<int, int> cellUV() const { return std::pair<int, int>(cellU(), cellV()); }
63  int cellX() const {
64  int N = (type() == HGCalFine) ? HGCalFineN : HGCalCoarseN;
65  return (3 * (cellV() - N) + 2);
66  }
67  int cellY() const {
68  int N = (type() == HGCalFine) ? HGCalFineN : HGCalCoarseN;
69  return (2 * cellU() - (N + cellV()));
70  }
71  std::pair<int, int> cellXY() const { return std::pair<int, int>(cellX(), cellY()); }
72 
74  int waferUAbs() const { return (id_ >> kHGCalWaferUOffset) & kHGCalWaferUMask; }
75  int waferVAbs() const { return (id_ >> kHGCalWaferVOffset) & kHGCalWaferVMask; }
76  int waferU() const { return (((id_ >> kHGCalWaferUSignOffset) & kHGCalWaferUSignMask) ? -waferUAbs() : waferUAbs()); }
77  int waferV() const { return (((id_ >> kHGCalWaferVSignOffset) & kHGCalWaferVSignMask) ? -waferVAbs() : waferVAbs()); }
78  std::pair<int, int> waferUV() const { return std::pair<int, int>(waferU(), waferV()); }
79  int waferX() const { return (-2 * waferU() + waferV()); }
80  int waferY() const { return (2 * waferV()); }
81  std::pair<int, int> waferXY() const { return std::pair<int, int>(waferX(), waferY()); }
82 
83  // get trigger cell u,v
84  int triggerCellU() const {
85  int N = (type() == HGCalFine) ? HGCalFineN : HGCalCoarseN;
87  return (cellU() >= N && cellV() >= N)
88  ? cellU() / NT
89  : ((cellU() < N && cellU() <= cellV()) ? cellU() / NT : (1 + (cellU() - (cellV() % NT + 1)) / NT));
90  }
91  int triggerCellV() const {
92  int N = (type() == HGCalFine) ? HGCalFineN : HGCalCoarseN;
94  return (cellU() >= N && cellV() >= N)
95  ? cellV() / NT
96  : ((cellU() < N && cellU() <= cellV()) ? ((cellV() - cellU()) / NT + cellU() / NT) : cellV() / NT);
97  }
98  std::pair<int, int> triggerCellUV() const { return std::pair<int, int>(triggerCellU(), triggerCellV()); }
99 
101  bool isEE() const { return (det() == HGCalEE); }
102  bool isHE() const { return (det() == HGCalHSi); }
103  bool isForward() const { return true; }
104 
106 
107 public:
108  static const int kHGCalCellUOffset = 0;
109  static const int kHGCalCellUMask = 0x1F;
110  static const int kHGCalCellVOffset = 5;
111  static const int kHGCalCellVMask = 0x1F;
112  static const int kHGCalWaferUOffset = 10;
113  static const int kHGCalWaferUMask = 0xF;
114  static const int kHGCalWaferUSignOffset = 14;
115  static const int kHGCalWaferUSignMask = 0x1;
116  static const int kHGCalWaferVOffset = 15;
117  static const int kHGCalWaferVMask = 0xF;
118  static const int kHGCalWaferVSignOffset = 19;
119  static const int kHGCalWaferVSignMask = 0x1;
120  static const int kHGCalLayerOffset = 20;
121  static const int kHGCalLayerMask = 0x1F;
122  static const int kHGCalZsideOffset = 25;
123  static const int kHGCalZsideMask = 0x1;
124  static const int kHGCalTypeOffset = 26;
125  static const int kHGCalTypeMask = 0x3;
126 };
127 
128 std::ostream& operator<<(std::ostream&, const HGCSiliconDetId& id);
129 
130 #endif
static const int kHGCalWaferVOffset
std::pair< int, int > cellUV() const
bool isForward() const
static const int kHGCalTypeMask
int waferX() const
int cellX() const
int waferU() const
int cellV() const
static const int kHGCalWaferVSignOffset
std::pair< int, int > waferUV() const
static const int kHGCalCellUMask
int zside() const
get the z-side of the cell (1/-1)
static const int HGCalFineN
std::ostream & operator<<(std::ostream &out, const ALILine &li)
Definition: ALILine.cc:167
static const int kHGCalCellVOffset
HGCSiliconDetId moduleId() const
std::pair< int, int > cellXY() const
int cellU() const
get the cell #&#39;s in u,v or in x,y
static const int kHGCalTypeOffset
DetId::Detector subdet() const
get the subdetector
HGCSiliconDetId geometryCell() const
static const int kHGCalCellVMask
static const int kHGCalZsideOffset
int type() const
get the type
std::pair< int, int > waferXY() const
int layer() const
get the layer #
static const int HGCalCoarseTrigger
static const HGCSiliconDetId Undefined
int waferV() const
bool isHE() const
int waferVAbs() const
std::pair< int, int > triggerCellUV() const
static const int HGCalCoarseN
static const int kHGCalLayerMask
bool isEE() const
consistency check : no bits left =&gt; no overhead
Definition: DetId.h:17
static const int kHGCalWaferVMask
HGCSiliconDetId & operator=(const DetId &id)
static const int kHGCalWaferUSignMask
static const int kHGCalWaferUMask
#define N
Definition: blowfish.cc:9
static const int kHGCalZsideMask
static const int kHGCalWaferVSignMask
Detector
Definition: DetId.h:24
uint32_t id_
Definition: DetId.h:69
static const int kHGCalCellUOffset
int waferUAbs() const
get the wafer #&#39;s in u,v or in x,y
int waferY() const
static const int kHGCalLayerOffset
int triggerCellU() const
int triggerCellV() const
static const int kHGCalWaferUOffset
static const int HGCalFineTrigger
static const int kHGCalWaferUSignOffset
int cellY() const
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:46