CMS 3D CMS Logo

HGCalCellOffset.cc
Go to the documentation of this file.
4 #include <vector>
5 #include <iostream>
6 
7 //#define EDM_ML_DEBUG
8 
10  double waferSize, int32_t nFine, int32_t nCoarse, double guardRingOffset_, double mouseBiteCut_) {
11  ncell_[0] = nFine;
12  ncell_[1] = nCoarse;
13  hgcalcell_ = std::make_unique<HGCalCell>(waferSize, nFine, nCoarse);
14  for (int k = 0; k < 2; ++k) { // k refers to type of wafer fine or coarse
15  cellX_[k] = waferSize / (3 * ncell_[k]);
16  cellY_[k] = sqrt3By2_ * cellX_[k];
17  // For formulas used please refer to https://indico.cern.ch/event/1297259/contributions/5455745/attachments/2667954/4722855/Cell_centroid.pdf
18  for (int j = 0; j < 4; ++j) { // j refers to type of cell : corner, truncated, extended
19  if (j == HGCalCell::fullCell) {
20  for (int i = 0; i < 6; ++i) {
21  offsetX[k][j][i] = 0.0;
22  offsetY[k][j][i] = 0.0;
23  cellArea[k][j] = 3 * sqrt3By2_ * cellY_[k];
24  }
25  } else if (j == HGCalCell::cornerCell) { // Offset for corner cells
26  double totalArea = (11.0 * sqrt3_ / 8.0) * std::pow(cellY_[k], 2); // Area of cell without any dead zone
27  double cutArea1 = (cellY_[k] * sqrt3_ * guardRingOffset_) - (std::pow(guardRingOffset_, 2) / (2 * sqrt3_));
28  double cutArea2 = (cellY_[k] * sqrt3By2_ * guardRingOffset_) - (std::pow(guardRingOffset_, 2) / (2 * sqrt3_));
29  double cutArea3 = sqrt3_ * std::pow((mouseBiteCut_ - (guardRingOffset_ / sqrt3By2_)), 2);
30  double x1_0 =
31  ((1.5 * cellY_[k]) - (cellY_[k] * 0.5 * guardRingOffset_) + (std::pow(guardRingOffset_, 2) / 18)) /
32  ((cellY_[k] * sqrt3_) - (guardRingOffset_ / (2 * sqrt3_)));
33  double y1_0 = ((cellY_[k] * -sqrt3By2_ * guardRingOffset_) + (std::pow(guardRingOffset_, 2) / (3 * sqrt3_))) /
34  ((cellY_[k] * sqrt3_) - (guardRingOffset_ / (2 * sqrt3_)));
35  double x2_0 =
36  ((0.375 * cellY_[k]) - (cellY_[k] * 0.25 * guardRingOffset_) + (std::pow(guardRingOffset_, 2) / 18)) /
37  ((cellY_[k] * sqrt3By2_) - (guardRingOffset_ / (2 * sqrt3_)));
38  double y2_0 =
39  ((cellY_[k] * -0.25 * sqrt3_ * guardRingOffset_) + (std::pow(guardRingOffset_, 2) / (3 * sqrt3_))) /
40  ((cellY_[k] * sqrt3By2_) - (guardRingOffset_ / (2 * sqrt3_)));
41  double x1 = (sqrt3By2_ * x1_0) - (0.5 * y1_0) - cellY_[k];
42  double y1 = (sqrt3By2_ * y1_0) + (0.5 * x1_0);
43  double x2 = (-sqrt3By2_ * x2_0) - (0.5 * y2_0) + (cellY_[k] * 1.25);
44  double y2 = (sqrt3By2_ * y2_0) + (0.5 * x2_0) + (cellY_[k] * 0.25 * sqrt3_);
45  double x3 = 0.5 * cellY_[k];
46  double y3 = (cellY_[k] * sqrt3By2_) - (2 * mouseBiteCut_ / 3) - (guardRingOffset_ / (3 * sqrt3By2_));
47  cellArea[k][j] = totalArea - cutArea1 - cutArea2 - cutArea3;
48  // Magnitude of offset of bottom corner cell in forward wafer
49  double xMag = ((35.0 * cellY_[k] / 132.0) * totalArea - (cutArea1 * x1) - (cutArea2 * x2) - (cutArea3 * x3)) /
50  (cellArea[k][j]);
51  double yMag =
52  ((5.0 * cellY_[k] / (44.0 * sqrt3_)) * totalArea - (cutArea1 * y1) - (cutArea2 * y2) - (cutArea3 * y3)) /
53  (cellArea[k][j]);
54  // (x, y) coordinates of offset for 6 corners of wafer starting from bottomCorner in clockwise direction
55  // offset_x = -1^(i%2)*(Offset_magnitude_X * cos(60*i) - Offset_magnitude_Y * sin(60*i)) i in (0-6)
56  // offset_x = (Offset_magnitude_Y * cos(60*i) + Offset_magnitude_X * sin(60*i)) i in (0-6)
57  std::array<double, 6> tempOffsetX = {{xMag,
58  -1.0 * (0.5 * xMag - sqrt3By2_ * yMag),
59  (-0.5 * xMag + sqrt3By2_ * yMag),
60  xMag,
61  (-0.5 * xMag - sqrt3By2_ * yMag),
62  -1.0 * (0.5 * xMag + sqrt3By2_ * yMag)}};
63  std::array<double, 6> tempOffsetY = {{yMag,
64  (0.5 * yMag + sqrt3By2_ * xMag),
65  (-0.5 * yMag - sqrt3By2_ * xMag),
66  -yMag,
67  (-0.5 * yMag + sqrt3By2_ * xMag),
68  (0.5 * yMag - sqrt3By2_ * xMag)}};
69  for (int i = 0; i < 6; ++i) {
70  offsetX[k][j][i] = tempOffsetX[i];
71  offsetY[k][j][i] = tempOffsetY[i];
72  }
73  } else if (j == HGCalCell::truncatedCell) { // Offset for truncated cells
74  double totalArea = (5.0 * sqrt3_ / 4.0) * std::pow(cellY_[k], 2); // Area of cell without any dead zone
75  double cutArea =
76  cellY_[k] * sqrt3_ * guardRingOffset_; // Area of inactive region form guardring and other effects
77  cellArea[k][j] = totalArea - cutArea;
78  double offMag = (((-2.0 / 15.0) * totalArea * cellY_[k]) - ((cellY_[k] - (0.5 * guardRingOffset_)) * cutArea)) /
79  (cellArea[k][j]); // Magnitude of offset
80  // (x, y) coordinates of offset for 6 sides of wafer starting from bottom left edge in clockwise direction
81  // offset_x = -Offset_magnitude * sin(30 + 60*i) i in (0-6)
82  // offset_y = -Offset_magnitude * cos(30 + 60*i) i in (0-6)
83  std::array<double, 6> tempOffsetX = {
84  {-0.5 * offMag, -offMag, -0.5 * offMag, 0.5 * offMag, offMag, 0.5 * offMag}};
85  std::array<double, 6> tempOffsetY = {
86  {-sqrt3By2_ * offMag, 0.0, sqrt3By2_ * offMag, sqrt3By2_ * offMag, 0.0, -sqrt3By2_ * offMag}};
87  for (int i = 0; i < 6; ++i) {
88  offsetX[k][j][i] = tempOffsetX[i];
89  offsetY[k][j][i] = tempOffsetY[i];
90  }
91  } else if (j == HGCalCell::extendedCell) { //Offset for extended cells
92  double totalArea = (7.0 * sqrt3_ / 4.0) * std::pow(cellY_[k], 2); // Area of cell without any dead zone
93  double cutArea =
94  cellY_[k] * sqrt3_ * guardRingOffset_; // Area of inactive region form guardring and other effects
95  cellArea[k][j] = totalArea - cutArea;
96  double offMag = // Magnitude of offset
97  (((5.0 / 42.0) * totalArea * cellY_[k]) - ((cellY_[k] - (0.5 * guardRingOffset_))) * (cutArea)) /
98  (cellArea[k][j]);
99  // (x, y) coordinates of offset for 6 sides of wafer starting from bottom left edge in clockwise direction
100  // offset_x = -Offset_magnitude * sin(30 + 60*i) i in (0-6)
101  // offset_y = -Offset_magnitude * cos(30 + 60*i) i in (0-6)
102  std::array<double, 6> tempOffsetX = {
103  {-0.5 * offMag, -offMag, -0.5 * offMag, 0.5 * offMag, offMag, 0.5 * offMag}};
104  std::array<double, 6> tempOffsetY = {
105  {-sqrt3By2_ * offMag, 0.0, sqrt3By2_ * offMag, sqrt3By2_ * offMag, 0.0, -sqrt3By2_ * offMag}};
106  for (int i = 0; i < 6; ++i) {
107  offsetX[k][j][i] = tempOffsetX[i];
108  offsetY[k][j][i] = tempOffsetY[i];
109  }
110  }
111  }
112  }
113 #ifdef EDM_ML_DEBUG
114  edm::LogVerbatim("HGCalGeom") << "HGCalCellOffset initialized with waferSize " << waferSize << " number of cells "
115  << nFine << ":" << nCoarse << " Guardring offset " << guardRingOffset_ << " Mousebite "
116  << mouseBiteCut_;
117 #endif
118 }
119 
120 std::pair<double, double> HGCalCellOffset::cellOffsetUV2XY1(int32_t u, int32_t v, int32_t placementIndex, int32_t type) {
121  if (type != 0)
122  type = 1;
123  double x_off(0), y_off(0);
124  std::pair<int, int> cell = hgcalcell_->cellType(u, v, ncell_[type], placementIndex);
125  int cellPos = cell.first;
126  int cellType = cell.second;
127  if (cellType == HGCalCell::truncatedCell || cellType == HGCalCell::extendedCell) {
128  x_off = offsetX[type][cellType][cellPos - HGCalCell::bottomLeftEdge];
129  y_off = offsetY[type][cellType][cellPos - HGCalCell::bottomLeftEdge];
130  } else if (cellType == HGCalCell::cornerCell) {
131  // The offset fo corner cells, is flipped along y-axis for 60 degree rotation of wafer
132  // and from forward to backward wafers
133  if (((placementIndex >= HGCalCell::cellPlacementExtra) && (placementIndex % 2 == 0)) ||
134  ((placementIndex < HGCalCell::cellPlacementExtra) && (placementIndex % 2 == 1))) {
135  cellPos = HGCalCell::bottomCorner + (6 + HGCalCell::bottomCorner - cellPos) % 6;
136  }
137  x_off = offsetX[type][cellType][cellPos - HGCalCell::bottomCorner];
138  y_off = offsetY[type][cellType][cellPos - HGCalCell::bottomCorner];
139  if (((placementIndex >= HGCalCell::cellPlacementExtra) && (placementIndex % 2 == 0)) ||
140  ((placementIndex < HGCalCell::cellPlacementExtra) && (placementIndex % 2 == 1))) {
141  x_off = -1 * x_off;
142  }
143  }
144 #ifdef EDM_ML_DEBUG
145  edm::LogVerbatim("HGCalGeom") << "HGCalCellOffset in wafer with placement index " << placementIndex << " type "
146  << type << " for cell u " << u << " v " << v << " Offset x:y " << x_off << ":" << y_off;
147 #endif
148  return std::make_pair(x_off, y_off);
149 }
150 
151 double HGCalCellOffset::cellAreaUV(int32_t u, int32_t v, int32_t placementIndex, int32_t type, bool reco) {
152  if (type != 0)
153  type = 1;
154  double area(0);
155  std::pair<int, int> cell = hgcalcell_->cellType(u, v, ncell_[type], placementIndex);
156  int cellType = cell.second;
157  area = reco ? cellArea[type][cellType] : HGCalParameters::k_ScaleToDDD2 * cellArea[type][cellType];
158 #ifdef EDM_ML_DEBUG
159  edm::LogVerbatim("HGCalGeom") << "HGCalCellOffset in wafer with placement index " << placementIndex << " type "
160  << type << " for cell u " << u << " v " << v << " Area " << area;
161 #endif
162  return area;
163 }
Log< level::Info, true > LogVerbatim
static constexpr int32_t fullCell
Definition: HGCalCell.h:28
double cellArea[2][4]
int32_t ncell_[2]
std::array< std::array< std::array< double, 6 >, 4 >, 2 > offsetY
std::pair< double, double > cellOffsetUV2XY1(int32_t u, int32_t v, int32_t placementIndex, int32_t type)
constexpr int pow(int x)
Definition: conifer.h:24
const double sqrt3By2_
double cellAreaUV(int32_t u, int32_t v, int32_t placementIndex, int32_t type, bool reco)
const double sqrt3_
static constexpr int32_t bottomLeftEdge
Definition: HGCalCell.h:35
std::array< std::array< std::array< double, 6 >, 4 >, 2 > offsetX
static constexpr int32_t truncatedCell
Definition: HGCalCell.h:30
static constexpr int32_t extendedCell
Definition: HGCalCell.h:31
std::unique_ptr< HGCalCell > hgcalcell_
fixed size matrix
HGCalCellOffset(double waferSize, int32_t nFine, int32_t nCoarse, double guardRingOffset_, double mouseBiteCut_)
static constexpr int32_t bottomCorner
Definition: HGCalCell.h:41
static constexpr int32_t cornerCell
Definition: HGCalCell.h:29
static constexpr double k_ScaleToDDD2
static constexpr int32_t cellPlacementExtra
Definition: HGCalCell.h:24