CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
HGCalCellOffset Class Reference

#include <HGCalCellOffset.h>

Public Member Functions

double cellAreaUV (int32_t u, int32_t v, int32_t placementIndex, int32_t type, bool reco)
 
double cellAreaUV (int32_t u, int32_t v, int32_t placementIndex, int32_t type, int32_t partialType, bool reco)
 
std::pair< double, double > cellOffsetUV2XY1 (int32_t u, int32_t v, int32_t placementIndex, int32_t type)
 
std::pair< double, double > cellOffsetUV2XY1 (int32_t u, int32_t v, int32_t placementIndex, int32_t type, int32_t partialType)
 
 HGCalCellOffset (double waferSize, int32_t nFine, int32_t nCoarse, double guardRingOffset_, double mouseBiteCut_, double sizeOffset_)
 

Private Attributes

double cellArea [2][6]
 
double cellAreaPartial [2][11]
 
double cellX_ [2]
 
double cellY_ [2]
 
double fullArea [2]
 
std::unique_ptr< HGCalCellhgcalcell_
 
int32_t ncell_ [2]
 
std::array< std::array< std::array< double, 6 >, 11 >, 2 > offsetPartialX
 
std::array< std::array< std::array< double, 6 >, 11 >, 2 > offsetPartialY
 
std::array< std::array< std::array< double, 6 >, 6 >, 2 > offsetX
 
std::array< std::array< std::array< double, 6 >, 6 >, 2 > offsetY
 
const double sqrt3_ = std::sqrt(3.0)
 
const double sqrt3By2_ = (0.5 * sqrt3_)
 

Detailed Description

Definition at line 10 of file HGCalCellOffset.h.

Constructor & Destructor Documentation

◆ HGCalCellOffset()

HGCalCellOffset::HGCalCellOffset ( double  waferSize,
int32_t  nFine,
int32_t  nCoarse,
double  guardRingOffset_,
double  mouseBiteCut_,
double  sizeOffset_ 
)

Definition at line 9 of file HGCalCellOffset.cc.

References cellArea, cellAreaPartial, cellX_, cellY_, HGCalCell::cornerCell, HGCalCell::extendedCell, HGCalCell::extendedMBCell, HGCalCell::fullCell, h, data-class-funcs::H, HGCalCell::halfCell, HGCalCell::halfExtCell, HGCalCell::halfTrunCell, HGCalCell::HDPartial0920Cell, HGCalCell::HDPartial1021Cell, hgcalcell_, mps_fire::i, dqmiolumiharvest::j, isotrackApplyRegressor::k, HGCalCell::LDPartial0007Cell, HGCalCell::LDPartial0209Cell, HGCalCell::LDPartial0714Cell, HGCalCell::LDPartial0815Cell, HGCalCell::LDPartial1415Cell, HGCalCell::LDPartial1515Cell, ncell_, offsetPartialX, offsetPartialY, offsetX, offsetY, HGCalCell::partiaclWaferCellsOffset, funct::pow(), sqrt3_, sqrt3By2_, HGCalCell::truncatedCell, HGCalCell::truncatedMBCell, testProducerWithPsetDescEmpty_cfi::x1, testProducerWithPsetDescEmpty_cfi::x2, testProducerWithPsetDescEmpty_cfi::y1, and testProducerWithPsetDescEmpty_cfi::y2.

14  {
15  ncell_[0] = nFine;
16  ncell_[1] = nCoarse;
17  hgcalcell_ = std::make_unique<HGCalCell>(waferSize, nFine, nCoarse);
18  double guardRingSizeOffset_ = guardRingOffset_ + sizeOffset_;
19  for (int k = 0; k < 2; ++k) { // k refers to type of wafer fine or coarse
20  cellX_[k] = waferSize / (3 * ncell_[k]);
21  cellY_[k] = sqrt3By2_ * cellX_[k];
22  // For formulas used please refer to https://indico.cern.ch/event/1297259/contributions/5455745/attachments/2667954/4722855/Cell_centroid.pdf
23  for (int j = 0; j < 6; ++j) { // j refers to type of cell : corner, truncated, extended, truncatedMB, extendedMB
24  if (j == HGCalCell::fullCell) {
25  for (int i = 0; i < 6; ++i) {
26  offsetX[k][j][i] = 0.0;
27  offsetY[k][j][i] = 0.0;
28  cellArea[k][j] = 3 * sqrt3By2_ * std::pow(cellX_[k], 2);
29  }
30  } else if (j == HGCalCell::cornerCell) { // Offset for corner cells
31  if (k == 0) {
32  double h = (mouseBiteCut_ - sqrt3By2_ * cellX_[k]);
33  double H = mouseBiteCut_ + guardRingOffset_ - (1 / sqrt3By2_ * guardRingOffset_);
34  double h1 = H - (sqrt3_ / 4 * cellX_[k]) + (guardRingSizeOffset_ / (2 * sqrt3_));
35  double h2 = H - (sqrt3_ / 2 * cellX_[k]) + (guardRingSizeOffset_ / (2 * sqrt3_));
36  double totalArea = 11.0 * sqrt3_ * std::pow(cellX_[k], 2) / 8.0;
37  double cutArea1 =
38  (sqrt3By2_ * cellX_[k] * guardRingSizeOffset_) - (0.5 / sqrt3_ * std::pow(guardRingSizeOffset_, 2));
39  double cutArea2 =
40  (sqrt3_ * cellX_[k] * guardRingSizeOffset_) - (0.5 / sqrt3_ * std::pow(guardRingSizeOffset_, 2));
41  double A1 = 2.0 * cellX_[k] * h - std::pow(h, 2) / (sqrt3_);
42  double A2 = sqrt3By2_ * cellX_[k] * cellX_[k];
43  double A3 = sqrt3By2_ * cellX_[k] * cellX_[k] / 4.0;
44  double cutArea3 =
45  sqrt3_ * std::pow(H, 2) - (1 / sqrt3By2_ * std::pow(h1, 2)) - (1 / sqrt3By2_ * std::pow(h2, 2));
46 #ifdef EDM_ML_DEBUG
47  edm::LogVerbatim("HGCalGeomX") << "h1 " << h1 << " h2 " << h2 << " H " << H << " cutarea1 " << cutArea1
48  << " cutarea2 " << cutArea2 << " cutarea3 " << cutArea3 << " " << cellX_[k]
49  << " " << guardRingSizeOffset_;
50 #endif
51  double x3_1 = -(((2.0 * std::pow(h, 3)) / (3.0 * sqrt3_) - cellX_[k] * std::pow(h, 2)) / A1);
52  double y3_1 = 0;
53  double x3_2 = -(sqrt3By2_ * cellX_[k] / 3);
54  double y3_2 = cellX_[k] / 6.0;
55  double x3_3 = -(cellX_[k] * sqrt3_ / 4.0);
56  double y3_3 = cellX_[k] * 11.0 / 12.0;
57 
58  double x1 = -(3.0 * sqrt3_ * cellX_[k] / 8.0 - sqrt3_ * guardRingOffset_ / 4.0);
59  double y1 = 5.0 * cellX_[k] / 12.0 - guardRingOffset_ / 4.0;
60  double x2 = -((0.5 * cellX_[k] - 0.5 * guardRingOffset_) * sqrt3By2_);
61  double y2 = -(0.5 * cellX_[k] - 0.5 * guardRingOffset_) * 0.5;
62  double x3 = (A1 * x3_1 + A2 * x3_2 + A3 * x3_3) / cutArea3;
63  double y3 = (A1 * y3_1 + A2 * y3_2 + A3 * y3_3) / cutArea3;
64  cellArea[k][j] = totalArea - cutArea1 - cutArea2 - cutArea3;
65  double xMag1 =
66  ((5.0 * sqrt3_ * cellX_[k] / 132.0) * totalArea - (cutArea1 * x1) - (cutArea2 * x2) - (cutArea3 * x3)) /
67  (cellArea[k][j]);
68  double yMag1 =
69  ((19.0 * cellX_[k] / 132.0) * totalArea - (cutArea1 * y1) - (cutArea2 * y2) - (cutArea3 * y3)) /
70  (cellArea[k][j]);
71 
72  double xMag = 0.5 * xMag1 + sqrt3By2_ * yMag1;
73  double yMag = sqrt3By2_ * xMag1 - 0.5 * yMag1;
74 
75  std::array<double, 6> tempOffsetX = {{(sqrt3By2_ * xMag - 0.5 * yMag),
76  yMag,
77  yMag,
78  (sqrt3By2_ * xMag - 0.5 * yMag),
79  (-sqrt3By2_ * xMag - 0.5 * yMag),
80  (-sqrt3By2_ * xMag - 0.5 * yMag)}};
81  std::array<double, 6> tempOffsetY = {{(0.5 * xMag + sqrt3By2_ * yMag),
82  xMag,
83  -xMag,
84  (-0.5 * xMag - sqrt3By2_ * yMag),
85  (0.5 * xMag - sqrt3By2_ * yMag),
86  (-0.5 * xMag + sqrt3By2_ * yMag)}};
87 
88  for (int i = 0; i < 6; ++i) {
89  offsetX[k][j][i] = tempOffsetX[i];
90  offsetY[k][j][i] = tempOffsetY[i];
91  }
92  } else if (k == 1) {
93  double h = (mouseBiteCut_ - guardRingOffset_) / sqrt3By2_ - cellX_[k] / 2;
94  double H = mouseBiteCut_ + guardRingOffset_ - (1 / sqrt3By2_ * guardRingOffset_);
95  double h1 = H - (sqrt3_ / 4 * cellX_[k]) + (guardRingSizeOffset_ / (2 * sqrt3_));
96  double h2 = H - (sqrt3_ / 2 * cellX_[k]) + (guardRingSizeOffset_ / (2 * sqrt3_));
97  double totalArea = 11.0 * sqrt3_ * std::pow(cellX_[k], 2) / 8.0;
98  double cutArea1 =
99  (sqrt3By2_ * cellX_[k] * guardRingSizeOffset_) - (0.5 / sqrt3_ * std::pow(guardRingSizeOffset_, 2));
100  double cutArea2 =
101  (sqrt3_ * cellX_[k] * guardRingSizeOffset_) - (0.5 / sqrt3_ * std::pow(guardRingSizeOffset_, 2));
102  double cutArea3 =
103  sqrt3_ * std::pow(H, 2) - (1 / sqrt3By2_ * std::pow(h1, 2)) - (1 / sqrt3By2_ * std::pow(h2, 2));
104  //double cutArea3 = sqrt3_ * std::pow((mouseBiteCut_ - guardRingOffset_), 2) - sqrt3By2_ * std::pow(h, 2);
105 
106  double x2_0 = (0.375 * cellX_[k] * cellX_[k] - (0.25 * cellX_[k] * guardRingOffset_) +
107  std::pow(guardRingOffset_, 2) / 18) /
108  (sqrt3By2_ * cellX_[k] + guardRingOffset_ / (2 * sqrt3_));
109  double y2_0 = (sqrt3_ * cellX_[k] * guardRingOffset_ / 4 + std::pow(guardRingOffset_, 2) / (6 * sqrt3_)) /
110  (sqrt3By2_ * cellX_[k] + guardRingOffset_ / (2 * sqrt3_));
111  double x3_1 = -(cellX_[k] - guardRingOffset_ - 2 * (mouseBiteCut_ - guardRingOffset_) / 3) * sqrt3By2_;
112  double y3_1 = 0.5 * (cellX_[k] - guardRingOffset_ - 2 * (mouseBiteCut_ - guardRingOffset_) / 3);
113  double x3_2 = -((3 * cellX_[k] / 2 - h / 3) * sqrt3By2_ + sqrt3_ * h / 6);
114  double y3_2 = -(cellX_[k] / 4 + 4 * h / 6);
115  double A1 = sqrt3_ * std::pow((mouseBiteCut_ - guardRingOffset_), 2);
116  double A2 = sqrt3By2_ * std::pow(h, 2);
117 
118  double x1 = 0;
119  double y1 = 0.5 * cellX_[k] - 0.5 * guardRingOffset_;
120  double x2 = -(1.5 * sqrt3By2_ * cellX_[k] - x2_0 * 0.5 - y2_0 * sqrt3By2_);
121  double y2 = -(0.25 * cellX_[k] - x2_0 * sqrt3By2_ + y2_0 / 2);
122  double x3 = (A1 * x3_1 - A2 * x3_2) / (A1 - A2);
123  double y3 = (A1 * y3_1 - A2 * y3_2) / (A1 - A2);
124  cellArea[k][j] = totalArea - cutArea1 - cutArea2 - cutArea3;
125  double xMag1 = ((0.0) * totalArea - (cutArea1 * x1) - (cutArea2 * x2) - (cutArea3 * x3)) / (cellArea[k][j]);
126  double yMag1 = ((-5 * cellX_[k] / 42) * totalArea - (cutArea1 * y1) - (cutArea2 * y2) - (cutArea3 * y3)) /
127  (cellArea[k][j]);
128 
129  double xMag = -0.5 * xMag1 - sqrt3By2_ * yMag1;
130  double yMag = sqrt3By2_ * xMag1 - 0.5 * yMag1;
131 
132  std::array<double, 6> tempOffsetX = {{(sqrt3By2_ * xMag - 0.5 * yMag),
133  yMag,
134  yMag,
135  (sqrt3By2_ * xMag - 0.5 * yMag),
136  (-sqrt3By2_ * xMag - 0.5 * yMag),
137  (-sqrt3By2_ * xMag - 0.5 * yMag)}};
138  std::array<double, 6> tempOffsetY = {{(0.5 * xMag + sqrt3By2_ * yMag),
139  xMag,
140  -xMag,
141  (-0.5 * xMag - sqrt3By2_ * yMag),
142  (0.5 * xMag - sqrt3By2_ * yMag),
143  (-0.5 * xMag + sqrt3By2_ * yMag)}};
144  for (int i = 0; i < 6; ++i) {
145  offsetX[k][j][i] = tempOffsetX[i];
146  offsetY[k][j][i] = tempOffsetY[i];
147  }
148  }
149  } else if (j == HGCalCell::truncatedCell) { // Offset for truncated cells
150  double totalArea = (5.0 * sqrt3_ / 4.0) * std::pow(cellX_[k], 2); // Area of cell without any dead zone
151  double cutArea =
152  cellX_[k] * sqrt3_ * guardRingSizeOffset_; // Area of inactive region form guardring and other effects
153  cellArea[k][j] = totalArea - cutArea;
154  double offMag = (((-2.0 / 15.0) * totalArea * cellX_[k]) - ((cellX_[k] - (0.5 * guardRingOffset_)) * cutArea)) /
155  (cellArea[k][j]); // Magnitude of offset
156  // (x, y) coordinates of offset for 6 sides of wafer starting from bottom left edge in clockwise direction
157  // offset_x = -Offset_magnitude * sin(30 + 60*i) i in (0-6)
158  // offset_y = -Offset_magnitude * cos(30 + 60*i) i in (0-6)
159  std::array<double, 6> tempOffsetX = {
160  {-0.5 * offMag, -offMag, -0.5 * offMag, 0.5 * offMag, offMag, 0.5 * offMag}};
161  std::array<double, 6> tempOffsetY = {
162  {-sqrt3By2_ * offMag, 0.0, sqrt3By2_ * offMag, sqrt3By2_ * offMag, 0.0, -sqrt3By2_ * offMag}};
163  for (int i = 0; i < 6; ++i) {
164  offsetX[k][j][i] = tempOffsetX[i];
165  offsetY[k][j][i] = tempOffsetY[i];
166  }
167  } else if (j == HGCalCell::extendedCell) { //Offset for extended cells
168  double totalArea = (7.0 * sqrt3_ / 4.0) * std::pow(cellX_[k], 2); // Area of cell without any dead zone
169  double cutArea =
170  cellX_[k] * sqrt3_ * guardRingSizeOffset_; // Area of inactive region form guardring and other effects
171  cellArea[k][j] = totalArea - cutArea;
172  double offMag = // Magnitude of offset
173  (((5.0 / 42.0) * totalArea * cellX_[k]) - ((cellX_[k] - (0.5 * guardRingOffset_))) * (cutArea)) /
174  (cellArea[k][j]);
175  // (x, y) coordinates of offset for 6 sides of wafer starting from bottom left edge in clockwise direction
176  // offset_x = -Offset_magnitude * sin(30 + 60*i) i in (0-6)
177  // offset_y = -Offset_magnitude * cos(30 + 60*i) i in (0-6)
178  std::array<double, 6> tempOffsetX = {
179  {-0.5 * offMag, -offMag, -0.5 * offMag, 0.5 * offMag, offMag, 0.5 * offMag}};
180  std::array<double, 6> tempOffsetY = {
181  {-sqrt3By2_ * offMag, 0.0, sqrt3By2_ * offMag, sqrt3By2_ * offMag, 0.0, -sqrt3By2_ * offMag}};
182  for (int i = 0; i < 6; ++i) {
183  offsetX[k][j][i] = tempOffsetX[i];
184  offsetY[k][j][i] = tempOffsetY[i];
185  }
186  } else if (j == HGCalCell::truncatedMBCell) {
187  double H = mouseBiteCut_ + guardRingOffset_ - (1 / sqrt3By2_ * guardRingOffset_);
188  double h = H - (sqrt3_ / 2 * cellX_[k]) + (guardRingSizeOffset_ / (2 * sqrt3_));
189  if (h > 0) {
190  double totalArea = 5.0 * sqrt3_ * std::pow(cellX_[k], 2) / 4.0;
191 
192  double cutArea1 = (sqrt3_ * cellX_[k] * guardRingSizeOffset_);
193  double cutArea2 = std::pow(h, 2) / sqrt3By2_;
194 
195  double x1 = -(0.5 * cellX_[k] - 0.5 * guardRingOffset_) * sqrt3By2_;
196  double y1 = -(0.5 * cellX_[k] - 0.5 * guardRingOffset_) * 0.5;
197  double x2 = -((sqrt3By2_ * cellX_[k]) - (2.0 * h) / 3.0);
198  double y2 = 0.5 * cellX_[k] - (2.0 * h) / (3.0 * sqrt3_);
199  cellArea[k][j] = totalArea - cutArea1 - cutArea2;
200 #ifdef EDM_ML_DEBUG
201  edm::LogVerbatim("HGCalGeomX") << "trunMB h " << h << " tot " << totalArea << " cutArea1 " << cutArea1
202  << " cutArea2 " << cutArea2;
203 #endif
204  double xMag1 =
205  ((sqrt3_ * cellX_[k] / 15.0) * totalArea - (cutArea1 * x1) - (cutArea2 * x2)) / (cellArea[k][j]);
206  double yMag1 = ((cellX_[k] / 15.0) * totalArea - (cutArea1 * y1) - (cutArea2 * y2)) / (cellArea[k][j]);
207  double xMag = -yMag1;
208  double yMag = -xMag1;
209 
210  std::array<double, 6> tempOffsetX = {{(sqrt3By2_ * xMag - 0.5 * yMag),
211  (-sqrt3By2_ * xMag - 0.5 * yMag),
212  (-sqrt3By2_ * xMag - 0.5 * yMag),
213  (sqrt3By2_ * xMag - 0.5 * yMag),
214  yMag,
215  yMag}};
216  std::array<double, 6> tempOffsetY = {{(0.5 * xMag + sqrt3By2_ * yMag),
217  (0.5 * xMag - sqrt3By2_ * yMag),
218  (-0.5 * xMag + sqrt3By2_ * yMag),
219  (-0.5 * xMag - sqrt3By2_ * yMag),
220  xMag,
221  -xMag}};
222  for (int i = 0; i < 6; ++i) {
223  offsetX[k][j][i] = tempOffsetX[i];
224  offsetY[k][j][i] = tempOffsetY[i];
225  }
226  } else {
228  std::array<double, 6> tempOffsetX = {{offsetX[k][HGCalCell::truncatedCell][0],
234  std::array<double, 6> tempOffsetY = {{offsetY[k][HGCalCell::truncatedCell][0],
240  for (int i = 0; i < 6; ++i) {
241  offsetX[k][j][i] = tempOffsetX[i];
242  offsetY[k][j][i] = tempOffsetY[i];
243  }
244  }
245  } else if (j == HGCalCell::extendedMBCell) {
246  double H = mouseBiteCut_ + guardRingOffset_ - (1 / sqrt3By2_ * guardRingOffset_);
247  double h = H - (sqrt3_ / 4 * cellX_[k]) + (guardRingSizeOffset_ / (2 * sqrt3_));
248 
249  double totalArea = 7.0 * sqrt3_ * std::pow(cellX_[k], 2) / 4.0;
250  double cutArea1 = (sqrt3_ * cellX_[k] * guardRingSizeOffset_);
251  double cutArea2 = std::pow(h, 2) / sqrt3By2_;
252 
253  double x1 = -(sqrt3By2_ * cellX_[k] - sqrt3By2_ * guardRingOffset_ / 2.0);
254  double y1 = (0.5 * cellX_[k] - 0.25 * guardRingOffset_);
255  double x2 = -(sqrt3By2_ * 1.5 * cellX_[k] - h / sqrt3_);
256  double y2 = -0.25 * cellX_[k] + h / 3.0;
257  cellArea[k][j] = totalArea - cutArea1 - cutArea2;
258 #ifdef EDM_ML_DEBUG
259  edm::LogVerbatim("HGCalGeomX") << H << "trunMB h " << h << " tot " << totalArea << " cutArea1 " << cutArea1
260  << " cutArea2 " << cutArea2;
261 #endif
262  double xMag1 =
263  ((-10.0 * sqrt3_ * cellX_[k] / 168.0) * totalArea - (cutArea1 * x1) - (cutArea2 * x2)) / (cellArea[k][j]);
264  double yMag1 = ((10.0 * cellX_[k] / 168.0) * totalArea - (cutArea1 * y1) - (cutArea2 * y2)) / (cellArea[k][j]);
265 
266  double xMag = yMag1;
267  double yMag = -xMag1;
268 
269  std::array<double, 6> tempOffsetX = {{(sqrt3By2_ * xMag - 0.5 * yMag),
270  yMag,
271  yMag,
272  (sqrt3By2_ * xMag - 0.5 * yMag),
273  (-sqrt3By2_ * xMag - 0.5 * yMag),
274  (-sqrt3By2_ * xMag - 0.5 * yMag)}};
275  std::array<double, 6> tempOffsetY = {{(0.5 * xMag + sqrt3By2_ * yMag),
276  xMag,
277  -xMag,
278  (-0.5 * xMag - sqrt3By2_ * yMag),
279  (0.5 * xMag - sqrt3By2_ * yMag),
280  (-0.5 * xMag + sqrt3By2_ * yMag)}};
281  for (int i = 0; i < 6; ++i) {
282  offsetX[k][j][i] = tempOffsetX[i];
283  offsetY[k][j][i] = tempOffsetY[i];
284  }
285  }
286  }
288  ++j) { //For cells in partial wafers
289  if (j == (HGCalCell::halfCell)) {
290  double totalArea = (3.0 * sqrt3_ / 4.0) * std::pow(cellX_[k], 2);
291  double cutArea = cellX_[k] * 2.0 * guardRingOffset_ - std::pow(guardRingOffset_, 2) / sqrt3_;
292  cellAreaPartial[k][j - HGCalCell::partiaclWaferCellsOffset] = totalArea - cutArea;
293  double x1 = (-cellX_[k] * guardRingOffset_ + 2 * std::pow(guardRingOffset_, 2) / (3 * sqrt3_)) /
294  (2 * cellX_[k] - guardRingOffset_ / sqrt3_);
295  double y1 = 0;
296  double xMag = ((-2.0 * sqrt3_ * cellX_[k] / 9.0) * totalArea - (cutArea * x1)) /
298  double yMag = (0 * totalArea - (cutArea * y1)) / (cellAreaPartial[k][j - HGCalCell::partiaclWaferCellsOffset]);
299 
300  std::array<double, 6> tempOffsetX = {{(-sqrt3By2_ * xMag - 0.5 * yMag),
301  (-sqrt3By2_ * xMag + 0.5 * yMag),
302  yMag,
303  (sqrt3By2_ * xMag + 0.5 * yMag),
304  (sqrt3By2_ * xMag - 0.5 * yMag),
305  -yMag}};
306  std::array<double, 6> tempOffsetY = {{(0.5 * xMag - sqrt3By2_ * yMag),
307  (-sqrt3By2_ * yMag - 0.5 * xMag),
308  -xMag,
309  (-0.5 * xMag + sqrt3By2_ * yMag),
310  (0.5 * xMag + sqrt3By2_ * yMag),
311  xMag}};
312 
313  for (int i = 0; i < 6; ++i) {
316  }
317  } else if (j == (HGCalCell::halfTrunCell)) {
318  double totalArea = 5 * sqrt3_ * std::pow(cellX_[k], 2) / 8;
319  double cutArea1 = (sqrt3By2_ * cellX_[k] * guardRingSizeOffset_) - guardRingOffset_ * guardRingSizeOffset_;
320  double cutArea2 = (3 * cellX_[k] * guardRingOffset_) / 2 - std::pow(guardRingOffset_, 2) / (2 * sqrt3_);
321 
322  double x1 = -sqrt3_ * cellX_[k] / 4;
323  double y1 = (0.5 * cellX_[k] - 0.5 * guardRingOffset_);
324  double x2 = (-3 * cellX_[k] * guardRingOffset_ / 4 + std::pow(guardRingOffset_, 2) / (3 * sqrt3_)) /
325  (3 * cellX_[k] / 2 - guardRingOffset_ / (2 * sqrt3_));
326  double y2 = (-cellX_[k] * guardRingOffset_ / (2 * sqrt3_) + std::pow(guardRingOffset_, 2) / 18 -
327  3 * std::pow(cellX_[k], 2) / 8) /
328  (3 * cellX_[k] / 2 - guardRingOffset_ / (2 * sqrt3_));
329  cellAreaPartial[k][j - HGCalCell::partiaclWaferCellsOffset] = totalArea - cutArea1 - cutArea2;
330  double xMag1 = ((-7 * sqrt3_ * cellX_[k] / 30) * totalArea - (cutArea1 * x1) - (cutArea2 * x2)) /
332  double yMag = ((-2 * cellX_[k] / 15) * totalArea - (cutArea1 * y1) - (cutArea2 * y2)) /
334  double xMag = -xMag1;
335 
336  std::array<double, 6> tempOffsetX = {{(-sqrt3By2_ * xMag - 0.5 * yMag),
337  (-sqrt3By2_ * xMag + 0.5 * yMag),
338  yMag,
339  (sqrt3By2_ * xMag + 0.5 * yMag),
340  (sqrt3By2_ * xMag - 0.5 * yMag),
341  -yMag}};
342  std::array<double, 6> tempOffsetY = {{(0.5 * xMag - sqrt3By2_ * yMag),
343  (-sqrt3By2_ * yMag - 0.5 * xMag),
344  -xMag,
345  (-0.5 * xMag + sqrt3By2_ * yMag),
346  (0.5 * xMag + sqrt3By2_ * yMag),
347  xMag}};
348  for (int i = 0; i < 6; ++i) {
351  }
352  } else if (j == (HGCalCell::halfExtCell)) {
353  double totalArea = (7.0 * sqrt3_ / 8.0) * std::pow(cellX_[k], 2);
354  double cutArea1 = cellX_[k] * sqrt3By2_ * guardRingOffset_ - std::pow(guardRingOffset_, 2);
355  double cutArea2 = cellX_[k] * 2.0 * guardRingOffset_ - std::pow(guardRingOffset_, 2) / (2 * sqrt3_);
356  cellAreaPartial[k][j - HGCalCell::partiaclWaferCellsOffset] = totalArea - cutArea1 - cutArea2;
357 
358  double x1 = -sqrt3By2_ * cellX_[k] / 2;
359  double y1 = -(cellX_[k] - guardRingOffset_ / 2);
360  double x2 = (-cellX_[k] * guardRingOffset_ + std::pow(guardRingOffset_, 2) / (3 * sqrt3_)) /
361  (2 * cellX_[k] - guardRingOffset_ / (2 * sqrt3_));
362  double y2 = (-cellX_[k] * guardRingOffset_ / (2 * sqrt3_) + std::pow(guardRingOffset_, 2) / 18) /
363  (2 * cellX_[k] - guardRingOffset_ / (2 * sqrt3_));
364  double xMag = ((-5 * sqrt3_ * cellX_[k] / 21.0) * totalArea - (cutArea1 * x1) - (cutArea2 * x2)) /
366  double yMag = ((-5 * cellX_[k] / 42.0) * totalArea - (cutArea1 * y1) - (cutArea2 * y2)) /
368 
369  std::array<double, 6> tempOffsetX = {{(-sqrt3By2_ * xMag - 0.5 * yMag),
370  (-sqrt3By2_ * xMag + 0.5 * yMag),
371  yMag,
372  (sqrt3By2_ * xMag + 0.5 * yMag),
373  (sqrt3By2_ * xMag - 0.5 * yMag),
374  -yMag}};
375  std::array<double, 6> tempOffsetY = {{(0.5 * xMag - sqrt3By2_ * yMag),
376  (-sqrt3By2_ * yMag - 0.5 * xMag),
377  -xMag,
378  (-0.5 * xMag + sqrt3By2_ * yMag),
379  (0.5 * xMag + sqrt3By2_ * yMag),
380  xMag}};
381 
382  for (int i = 0; i < 6; ++i) {
385  }
386  } else if (j == (HGCalCell::LDPartial0714Cell)) {
387  if (k == 1) {
388  double totalArea = (9.0 * sqrt3_ / 4.0) * std::pow(cellX_[k], 2);
389  double cutArea1 =
390  (3 * cellX_[k] * sqrt3By2_ * guardRingOffset_) - (std::pow(guardRingOffset_, 2) / (2 * sqrt3_));
391  double cutArea2 = (3 * cellX_[k] * sqrt3By2_ * guardRingOffset_);
392  double cutArea3 = sqrt3_ * std::pow((mouseBiteCut_ - (guardRingOffset_ / sqrt3By2_)), 2) / 2;
393  double x1_0 = ((3.375 * cellX_[k] * cellX_[k]) - (cellX_[k] * 0.75 * guardRingOffset_) +
394  (std::pow(guardRingOffset_, 2) / 18)) /
395  ((3 * cellX_[k] * sqrt3By2_) - (guardRingOffset_ / (2 * sqrt3_)));
396  double y1_0 =
397  ((3 * cellX_[k] * sqrt3By2_ * guardRingOffset_ / 2) - (std::pow(guardRingOffset_, 2) / (6 * sqrt3_))) /
398  ((3 * cellX_[k] * sqrt3By2_) - (guardRingOffset_ / (2 * sqrt3_)));
399 
400  double x2_0 = (3 * sqrt3By2_ * cellX_[k] / 2);
401  double y2_0 = guardRingOffset_ / 2;
402 
403  double x1 = (cellX_[k] / 2 - guardRingOffset_) * sqrt3By2_ + x1_0 * 0.5 + y1_0 * sqrt3By2_;
404  double y1 = cellX_[k] + (cellX_[k] / 2 - guardRingOffset_) * 0.5 - x1_0 * sqrt3By2_ + y1_0 * 0.5;
405 
406  double x2 = x2_0 - sqrt3By2_ * cellX_[k];
407  double y2 = -(cellX_[k] - y2_0);
408 
409  double x3 = sqrt3_ * cellX_[k] - mouseBiteCut_ + (mouseBiteCut_ - (guardRingOffset_ / sqrt3By2_)) / 3;
410  double y3 = -(cellX_[k] - sqrt3_ * (mouseBiteCut_ - (guardRingOffset_ / sqrt3By2_)) / 3 - guardRingOffset_);
411 
412  cellAreaPartial[k][j - HGCalCell::partiaclWaferCellsOffset] = totalArea - cutArea1 - cutArea2 - cutArea3;
413  double xMag = ((sqrt3_ * cellX_[k] / 8) * totalArea - (cutArea1 * x1) - (cutArea2 * x2) - (cutArea3 * x3)) /
415  double yMag = ((-1 * cellX_[k] / 8) * totalArea - (cutArea1 * y1) - (cutArea2 * y2) - (cutArea3 * y3)) /
417 
418  std::array<double, 6> tempOffsetX = {{(-sqrt3By2_ * xMag - 0.5 * yMag),
419  (-sqrt3By2_ * xMag + 0.5 * yMag),
420  yMag,
421  (sqrt3By2_ * xMag + 0.5 * yMag),
422  (sqrt3By2_ * xMag - 0.5 * yMag),
423  -yMag}};
424  std::array<double, 6> tempOffsetY = {{(0.5 * xMag - sqrt3By2_ * yMag),
425  (-sqrt3By2_ * yMag - 0.5 * xMag),
426  -xMag,
427  (-0.5 * xMag + sqrt3By2_ * yMag),
428  (0.5 * xMag + sqrt3By2_ * yMag),
429  xMag}};
430 
431  for (int i = 0; i < 6; ++i) {
434  }
435  } else {
437  for (int i = 0; i < 6; ++i) {
440  }
441  }
442  } else if (j == (HGCalCell::LDPartial0209Cell)) {
443  if (k == 1) {
444  double totalArea = (23.0 * sqrt3_ / 8.0) * std::pow(cellX_[k], 2);
445  double cutArea1 =
446  (5 * cellX_[k] * sqrt3By2_ * guardRingOffset_) - (std::pow(guardRingOffset_, 2) / (2 * sqrt3_));
447  double cutArea2 = (4 * cellX_[k] * guardRingOffset_);
448  double cutArea3 = std::pow(mouseBiteCut_, 2) / sqrt3_;
449 
450  double x1_0 = (9.375 * cellX_[k] * cellX_[k] - (cellX_[k] * 1.25 * guardRingOffset_) +
451  (std::pow(guardRingOffset_, 2) / 18)) /
452  ((5 * cellX_[k] * sqrt3By2_) - (guardRingOffset_ / (2 * sqrt3_)));
453  double y1_0 =
454  ((5 * cellX_[k] * sqrt3By2_ * guardRingOffset_ / 2) - (std::pow(guardRingOffset_, 2) / (6 * sqrt3_))) /
455  ((5 * cellX_[k] * sqrt3By2_) - (guardRingOffset_ / (2 * sqrt3_)));
456 
457  double x1 = (1.5 * cellX_[k]) * sqrt3By2_ - x1_0 * 0.5 - y1_0 * sqrt3By2_;
458  double y1 = -0.25 * cellX_[k] + x1_0 * sqrt3By2_ - y1_0 * 0.5;
459  double x2 = -(sqrt3By2_ * cellX_[k] - 0.5 * guardRingOffset_);
460  double y2 = 1.5 * cellX_[k];
461  double x3 = -(sqrt3By2_ * cellX_[k] - mouseBiteCut_ / 3);
462  double y3 = 3.5 * cellX_[k] - (5 * mouseBiteCut_) / 3 * sqrt3_;
463 
464  cellAreaPartial[k][j - HGCalCell::partiaclWaferCellsOffset] = totalArea - cutArea1 - cutArea2 - cutArea3;
465  double xMag =
466  ((-9 * cellX_[k] / (sqrt3_ * 92)) * totalArea - (cutArea1 * x1) - (cutArea2 * x2) - (cutArea3 * x3)) /
468  double yMag =
469  ((199 * cellX_[k] / (sqrt3_ * 276)) * totalArea - (cutArea1 * y1) - (cutArea2 * y2) - (cutArea3 * y3)) /
471 
472  std::array<double, 6> tempOffsetX = {{(-sqrt3By2_ * xMag - 0.5 * yMag),
473  (-sqrt3By2_ * xMag + 0.5 * yMag),
474  yMag,
475  (sqrt3By2_ * xMag + 0.5 * yMag),
476  (sqrt3By2_ * xMag - 0.5 * yMag),
477  -yMag}};
478  std::array<double, 6> tempOffsetY = {{(0.5 * xMag - sqrt3By2_ * yMag),
479  (-sqrt3By2_ * yMag - 0.5 * xMag),
480  -xMag,
481  (-0.5 * xMag + sqrt3By2_ * yMag),
482  (0.5 * xMag + sqrt3By2_ * yMag),
483  xMag}};
484  for (int i = 0; i < 6; ++i) {
487  }
488  } else {
490  for (int i = 0; i < 6; ++i) {
493  }
494  }
495  } else if (j == (HGCalCell::LDPartial0007Cell)) {
496  if (k == 1) {
497  double totalArea = (5.0 * sqrt3_ / 4.0) * std::pow(cellX_[k], 2);
498  double cutArea1 = (cellX_[k] * guardRingOffset_);
499  double cutArea2 = (sqrt3_ * cellX_[k] * guardRingOffset_);
500  double h = cellX_[k] - (sqrt3By2_ * cellX_[k] - mouseBiteCut_) / sqrt3By2_;
501  double cutArea3 = sqrt3_ * std::pow(h, 2) / 2;
502 
503  double x1 = cellX_[k] * sqrt3By2_ - guardRingOffset_ / 2;
504  double y1 = 0;
505  double x2 = 0;
506  double y2 = 0.5 * cellX_[k] - guardRingOffset_ / 2;
507  double x3 = sqrt3By2_ * cellX_[k] - guardRingOffset_ - h / sqrt3_;
508  double y3 = 0.5 * cellX_[k] - guardRingOffset_ - h / 3;
509 
510  cellAreaPartial[k][j - HGCalCell::partiaclWaferCellsOffset] = totalArea - cutArea1 - cutArea2 - cutArea3;
511  double xMag = ((0.0) * totalArea - (cutArea1 * x1) - (cutArea2 * x2) - (cutArea3 * x3)) /
513  double yMag = ((-2 * cellX_[k] / 15) * totalArea - (cutArea1 * y1) - (cutArea2 * y2) - (cutArea3 * y3)) /
515 
516  std::array<double, 6> tempOffsetX = {{(-sqrt3By2_ * xMag - 0.5 * yMag),
517  (-sqrt3By2_ * xMag + 0.5 * yMag),
518  yMag,
519  (sqrt3By2_ * xMag + 0.5 * yMag),
520  (sqrt3By2_ * xMag - 0.5 * yMag),
521  -yMag}};
522  std::array<double, 6> tempOffsetY = {{(0.5 * xMag - sqrt3By2_ * yMag),
523  (-sqrt3By2_ * yMag - 0.5 * xMag),
524  -xMag,
525  (-0.5 * xMag + sqrt3By2_ * yMag),
526  (0.5 * xMag + sqrt3By2_ * yMag),
527  xMag}};
528  for (int i = 0; i < 6; ++i) {
531  }
532  } else {
534  for (int i = 0; i < 6; ++i) {
537  }
538  }
539  } else if (j == (HGCalCell::LDPartial0815Cell)) {
540  if (k == 1) {
541  double totalArea = sqrt3_ * std::pow(cellX_[k], 2);
542  double cutArea1 = (sqrt3_ * cellX_[k] * guardRingOffset_);
543  double cutArea2 = (sqrt3_ * cellX_[k] * guardRingOffset_) - std::pow(guardRingOffset_, 2) / (2 * sqrt3_);
544  double cutArea3 = sqrt3_ * std::pow((mouseBiteCut_ - guardRingOffset_ / sqrt3By2_), 2) / 2;
545 
546  double x2_0 = (1.5 * cellX_[k] * cellX_[k] - (0.5 * cellX_[k] * guardRingOffset_) +
547  std::pow(guardRingOffset_, 2) / 18) /
548  (sqrt3_ * cellX_[k] - guardRingOffset_ / (2 * sqrt3_));
549  double y2_0 = (sqrt3By2_ * cellX_[k] * guardRingOffset_ - std::pow(guardRingOffset_, 2) / (sqrt3_ * 3)) /
550  (sqrt3_ * cellX_[k] - guardRingOffset_ / (2 * sqrt3_));
551  double x1 = 0;
552  double y1 = 0.5 * cellX_[k] - guardRingOffset_ / 2;
553  double x2 = x2_0 * 0.5 - y2_0 * sqrt3By2_;
554  double y2 = -(cellX_[k] - (x2_0 * sqrt3By2_ + y2_0 * 0.5));
555  double x3 = sqrt3By2_ * cellX_[k] - mouseBiteCut_ + (mouseBiteCut_ - guardRingOffset_ / sqrt3By2_) / 3;
556  double y3 = cellX_[k] - (mouseBiteCut_ - guardRingOffset_ / sqrt3By2_) / sqrt3_ - guardRingOffset_;
557 
558  cellAreaPartial[k][j - HGCalCell::partiaclWaferCellsOffset] = totalArea - cutArea1 - cutArea2 - cutArea3;
559  double xMag = -((-sqrt3_ * cellX_[k] / 9) * totalArea - (cutArea1 * x1) - (cutArea2 * x2) - (cutArea3 * x3)) /
561  double yMag = ((-cellX_[k] / 15) * totalArea - (cutArea1 * y1) - (cutArea2 * y2) - (cutArea3 * y3)) /
563 
564  std::array<double, 6> tempOffsetX = {{(-sqrt3By2_ * xMag - 0.5 * yMag),
565  (-sqrt3By2_ * xMag + 0.5 * yMag),
566  yMag,
567  (sqrt3By2_ * xMag + 0.5 * yMag),
568  (sqrt3By2_ * xMag - 0.5 * yMag),
569  -yMag}};
570  std::array<double, 6> tempOffsetY = {{(0.5 * xMag - sqrt3By2_ * yMag),
571  (-sqrt3By2_ * yMag - 0.5 * xMag),
572  -xMag,
573  (-0.5 * xMag + sqrt3By2_ * yMag),
574  (0.5 * xMag + sqrt3By2_ * yMag),
575  xMag}};
576  for (int i = 0; i < 6; ++i) {
579  }
580  } else {
581  for (int i = 0; i < 6; ++i) {
585  }
586  }
587  } else if (j == (HGCalCell::LDPartial1415Cell)) {
588  if (k == 1) {
589  double totalArea = 7 * sqrt3_ * std::pow(cellX_[k], 2) / 4;
590  double cutArea1 = (3 * cellX_[k] * guardRingOffset_);
591  double cutArea2 = (2 * sqrt3_ * cellX_[k] * guardRingOffset_) - std::pow(guardRingOffset_, 2) * sqrt3By2_;
592  double cutArea3 = std::pow((mouseBiteCut_ - guardRingOffset_), 2) / sqrt3_;
593 
594  double x2_0 = (6 * cellX_[k] * cellX_[k] - std::pow(guardRingOffset_, 2)) /
595  (2 * sqrt3_ * cellX_[k] - guardRingOffset_ * sqrt3By2_);
596  double y2_0 = (sqrt3_ * cellX_[k] * guardRingOffset_ - std::pow(guardRingOffset_, 2) / sqrt3_) /
597  (2 * sqrt3_ * cellX_[k] - guardRingOffset_ * sqrt3By2_);
598  double x1 = -sqrt3By2_ * cellX_[k] + guardRingOffset_ / 2;
599  double y1 = -cellX_[k];
600  double x2 = sqrt3By2_ * cellX_[k] - x2_0 * 0.5 - y2_0 * sqrt3By2_;
601  double y2 = 0.5 * cellX_[k] - x2_0 * sqrt3By2_ + y2_0 * 0.5;
602  double h = (mouseBiteCut_ - guardRingOffset_) / sqrt3By2_;
603  double x3 = -(cellX_[k] - h / 3 - guardRingOffset_) * sqrt3By2_;
604  double y3 = 5 * h / 6 - 5 * cellX_[k] / 2;
605 
606  cellAreaPartial[k][j - HGCalCell::partiaclWaferCellsOffset] = totalArea - cutArea1 - cutArea2 - cutArea3;
607  double xMag =
608  ((-2 * cellX_[k] / (7 * sqrt3_)) * totalArea - (cutArea1 * x1) - (cutArea2 * x2) - (cutArea3 * x3)) /
610  double yMag = ((-cellX_[k] / 3) * totalArea - (cutArea1 * y1) - (cutArea2 * y2) - (cutArea3 * y3)) /
612 
613  std::array<double, 6> tempOffsetX = {{(-sqrt3By2_ * xMag - 0.5 * yMag),
614  (-sqrt3By2_ * xMag + 0.5 * yMag),
615  yMag,
616  (sqrt3By2_ * xMag + 0.5 * yMag),
617  (sqrt3By2_ * xMag - 0.5 * yMag),
618  -yMag}};
619  std::array<double, 6> tempOffsetY = {{(0.5 * xMag - sqrt3By2_ * yMag),
620  (-sqrt3By2_ * yMag - 0.5 * xMag),
621  -xMag,
622  (-0.5 * xMag + sqrt3By2_ * yMag),
623  (0.5 * xMag + sqrt3By2_ * yMag),
624  xMag}};
625  for (int i = 0; i < 6; ++i) {
628  }
629  } else {
630  for (int i = 0; i < 6; ++i) {
634  }
635  }
636  } else if (j == (HGCalCell::LDPartial1515Cell)) {
637  if (k == 1) {
638  double totalArea = 7 * sqrt3_ * std::pow(cellX_[k], 2) / 8;
639  double cutArea1 = (2 * cellX_[k] * guardRingOffset_);
640  double cutArea2 = (sqrt3By2_ * cellX_[k] * guardRingOffset_);
641  double cutArea3 = cellX_[k] * (mouseBiteCut_ - guardRingOffset_) - (sqrt3_ * cellX_[k] * cellX_[k] / 8);
642 
643  double x1 = -guardRingOffset_ / 2;
644  double y1 = 0;
645  double x2 = -(sqrt3By2_ * cellX_[k] / 2);
646  double y2 = -(cellX_[k] - 0.5 * guardRingOffset_);
647  double x3 = (cellX_[k] * cellX_[k] / 8 - sqrt3_ * cellX_[k] * (mouseBiteCut_ - guardRingOffset_) / 4) /
648  ((mouseBiteCut_ - guardRingOffset_) - sqrt3_ * cellX_[k] / 8);
649  double y3 =
650  (std::pow((mouseBiteCut_ - guardRingOffset_), 2) / sqrt3_ -
651  (1.25 * cellX_[k] * (mouseBiteCut_ - guardRingOffset_)) + 7 * sqrt3_ * cellX_[k] * cellX_[k] / 48) /
652  ((mouseBiteCut_ - guardRingOffset_) - sqrt3_ * cellX_[k] / 8);
653 
654  cellAreaPartial[k][j - HGCalCell::partiaclWaferCellsOffset] = totalArea - cutArea1 - cutArea2 - cutArea3;
655  double xMag = (-(cellX_[k] / (sqrt3_)) * totalArea - (cutArea1 * x1) - (cutArea2 * x2) - (cutArea3 * x3)) /
657  double yMag = ((-5 * cellX_[k] / 42) * totalArea - (cutArea1 * y1) - (cutArea2 * y2) - (cutArea3 * y3)) /
659 
660  std::array<double, 6> tempOffsetX = {{(-sqrt3By2_ * xMag - 0.5 * yMag),
661  (-sqrt3By2_ * xMag + 0.5 * yMag),
662  yMag,
663  (sqrt3By2_ * xMag + 0.5 * yMag),
664  (sqrt3By2_ * xMag - 0.5 * yMag),
665  -yMag}};
666  std::array<double, 6> tempOffsetY = {{(0.5 * xMag - sqrt3By2_ * yMag),
667  (-sqrt3By2_ * yMag - 0.5 * xMag),
668  -xMag,
669  (-0.5 * xMag + sqrt3By2_ * yMag),
670  (0.5 * xMag + sqrt3By2_ * yMag),
671  xMag}};
672  for (int i = 0; i < 6; ++i) {
675  }
676  } else {
677  for (int i = 0; i < 6; ++i) {
681  }
682  }
683  } else if (j == (HGCalCell::HDPartial0920Cell)) {
684  if (k == 0) {
685  double totalArea = 37 * sqrt3_ * std::pow(cellX_[k], 2) / 24;
686  double cutArea1 = (4 * cellX_[k] * guardRingOffset_) / sqrt3_;
687  double cutArea2 =
688  (7 * cellX_[k] * guardRingOffset_) / (2 * sqrt3_) - std::pow(guardRingOffset_, 2) / (2 * sqrt3_);
689 
690  double x1 = cellX_[k] / (2 * sqrt3_);
691  double y1 = -(0.5 * cellX_[k] - 0.5 * guardRingOffset_);
692  double x2_0 = ((2.041 * cellX_[k] * cellX_[k]) - (cellX_[k] * 0.583 * guardRingOffset_) +
693  (std::pow(guardRingOffset_, 2) / 18)) /
694  ((7 * cellX_[k] / (2 * sqrt3_)) - (guardRingOffset_ / (2 * sqrt3_)));
695  double y2_0 =
696  ((7 * cellX_[k] * guardRingOffset_ / (4 * sqrt3_)) - std::pow(guardRingOffset_, 2) / (6 * sqrt3_)) /
697  ((7 * cellX_[k] / (2 * sqrt3_)) - (guardRingOffset_ / (2 * sqrt3_)));
698 
699  double x2 = (0.5 * x2_0) - (sqrt3By2_ * y2_0) + (cellX_[k] * 0.5 * sqrt3By2_);
700  double y2 = -(0.5 * y2_0) - (sqrt3By2_ * x2_0) + (cellX_[k] * 1.25);
701  cellAreaPartial[k][j - HGCalCell::partiaclWaferCellsOffset] = totalArea - cutArea1 - cutArea2;
702  double xMag = ((25 * sqrt3_ * cellX_[k] / 148) * totalArea - (cutArea1 * x1) - (cutArea2 * x2)) /
704  double yMag = ((73 * cellX_[k] / 444) * totalArea - (cutArea1 * y1) - (cutArea2 * y2)) /
706 
707  std::array<double, 6> tempOffsetX = {{(-sqrt3By2_ * xMag - 0.5 * yMag),
708  (-sqrt3By2_ * xMag + 0.5 * yMag),
709  yMag,
710  (sqrt3By2_ * xMag + 0.5 * yMag),
711  (sqrt3By2_ * xMag - 0.5 * yMag),
712  -yMag}};
713  std::array<double, 6> tempOffsetY = {{(0.5 * xMag - sqrt3By2_ * yMag),
714  (-sqrt3By2_ * yMag - 0.5 * xMag),
715  -xMag,
716  (-0.5 * xMag + sqrt3By2_ * yMag),
717  (0.5 * xMag + sqrt3By2_ * yMag),
718  xMag}};
719  for (int i = 0; i < 6; ++i) {
722  }
723  } else {
724  for (int i = 0; i < 6; ++i) {
728  }
729  }
730  } else if (j == (HGCalCell::HDPartial1021Cell)) {
731  if (k == 0) {
732  double totalArea = 11 * sqrt3_ * std::pow(cellX_[k], 2) / 6;
733  double cutArea1 = (5 * cellX_[k] * guardRingOffset_) / (2 * sqrt3_);
734  double cutArea2 =
735  (5 * cellX_[k] * guardRingOffset_) / (2 * sqrt3_) - std::pow(guardRingOffset_, 2) / (2 * sqrt3_);
736 
737  double x1 = -cellX_[k] / (4 * sqrt3_);
738  double y1 = cellX_[k] - 0.5 * guardRingOffset_;
739  double x2_0 = ((1.041 * cellX_[k] * cellX_[k]) - (cellX_[k] * 0.416 * guardRingOffset_) +
740  (std::pow(guardRingOffset_, 2) / 18.0)) /
741  ((5.0 * cellX_[k] / (2.0 * sqrt3_)) - (guardRingOffset_ / (2.0 * sqrt3_)));
742  double y2_0 =
743  ((5.0 * cellX_[k] * guardRingOffset_ / (4.0 * sqrt3_)) - std::pow(guardRingOffset_, 2) / (6 * sqrt3_)) /
744  ((5.0 * cellX_[k] / (2.0 * sqrt3_)) - (guardRingOffset_ / (2.0 * sqrt3_)));
745 
746  double x2 = -(0.5 * x2_0) + (sqrt3By2_ * y2_0) + (cellX_[k] * 1.5 * sqrt3By2_);
747  double y2 = -(0.5 * y2_0) + (sqrt3By2_ * x2_0) - cellX_[k];
748  cellAreaPartial[k][j - HGCalCell::partiaclWaferCellsOffset] = totalArea - cutArea1 - cutArea2;
749  double xMag = ((47.0 * cellX_[k] / (528.0 * sqrt3_)) * totalArea - (cutArea1 * x1) - (cutArea2 * x2)) /
751  double yMag = ((47.0 * cellX_[k] / 528.0) * totalArea - (cutArea1 * y1) - (cutArea2 * y2)) /
753 
754  std::array<double, 6> tempOffsetX = {{(-sqrt3By2_ * xMag - 0.5 * yMag),
755  (-sqrt3By2_ * xMag + 0.5 * yMag),
756  yMag,
757  (sqrt3By2_ * xMag + 0.5 * yMag),
758  (sqrt3By2_ * xMag - 0.5 * yMag),
759  -yMag}};
760  std::array<double, 6> tempOffsetY = {{(0.5 * xMag - sqrt3By2_ * yMag),
761  (-sqrt3By2_ * yMag - 0.5 * xMag),
762  -xMag,
763  (-0.5 * xMag + sqrt3By2_ * yMag),
764  (0.5 * xMag + sqrt3By2_ * yMag),
765  xMag}};
766  for (int i = 0; i < 6; ++i) {
769  }
770  } else {
771  for (int i = 0; i < 6; ++i) {
775  }
776  }
777  }
778  }
779  }
780 
781 #ifdef EDM_ML_DEBUG
782  edm::LogVerbatim("HGCalGeom") << "HGCalCellOffset initialized with waferSize " << waferSize << " number of cells "
783  << nFine << ":" << nCoarse << " Guardring offset " << guardRingOffset_ << " Mousebite "
784  << mouseBiteCut_;
785 #endif
786 }
static constexpr int32_t LDPartial0007Cell
Definition: HGCalCell.h:43
Log< level::Info, true > LogVerbatim
static constexpr int32_t fullCell
Definition: HGCalCell.h:28
int32_t ncell_[2]
double cellArea[2][6]
std::array< std::array< std::array< double, 6 >, 11 >, 2 > offsetPartialX
static constexpr int32_t HDPartial1021Cell
Definition: HGCalCell.h:49
static constexpr int32_t extendedMBCell
Definition: HGCalCell.h:33
static constexpr int32_t HDPartial0920Cell
Definition: HGCalCell.h:48
double cellAreaPartial[2][11]
static constexpr int32_t partiaclWaferCellsOffset
Definition: HGCalCell.h:39
static constexpr int32_t truncatedMBCell
Definition: HGCalCell.h:32
const double sqrt3By2_
static constexpr int32_t LDPartial0209Cell
Definition: HGCalCell.h:42
static constexpr int32_t halfCell
Definition: HGCalCell.h:36
static constexpr int32_t LDPartial0714Cell
Definition: HGCalCell.h:41
static constexpr int32_t LDPartial1515Cell
Definition: HGCalCell.h:46
const double sqrt3_
std::array< std::array< std::array< double, 6 >, 6 >, 2 > offsetY
static constexpr int32_t truncatedCell
Definition: HGCalCell.h:30
static constexpr int32_t halfTrunCell
Definition: HGCalCell.h:37
static constexpr int32_t extendedCell
Definition: HGCalCell.h:31
std::unique_ptr< HGCalCell > hgcalcell_
static constexpr int32_t LDPartial0815Cell
Definition: HGCalCell.h:44
static constexpr int32_t LDPartial1415Cell
Definition: HGCalCell.h:45
static constexpr int32_t halfExtCell
Definition: HGCalCell.h:38
static constexpr int32_t cornerCell
Definition: HGCalCell.h:29
std::array< std::array< std::array< double, 6 >, 11 >, 2 > offsetPartialY
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
std::array< std::array< std::array< double, 6 >, 6 >, 2 > offsetX
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29

Member Function Documentation

◆ cellAreaUV() [1/2]

double HGCalCellOffset::cellAreaUV ( int32_t  u,
int32_t  v,
int32_t  placementIndex,
int32_t  type,
bool  reco 
)

Definition at line 873 of file HGCalCellOffset.cc.

References custom_jme_cff::area, cellArea, hgcalcell_, HGCalParameters::k_ScaleToDDD2, ncell_, and findQualityFiles::v.

Referenced by cellAreaUV().

873  {
874  if (type != 0)
875  type = 1;
876  double area(0);
877  std::pair<int, int> cell = hgcalcell_->cellType(u, v, ncell_[type], placementIndex);
878  int cellType = cell.second;
879  area = reco ? cellArea[type][cellType] : HGCalParameters::k_ScaleToDDD2 * cellArea[type][cellType];
880  return area;
881 }
int32_t ncell_[2]
double cellArea[2][6]
std::unique_ptr< HGCalCell > hgcalcell_
fixed size matrix
static constexpr double k_ScaleToDDD2

◆ cellAreaUV() [2/2]

double HGCalCellOffset::cellAreaUV ( int32_t  u,
int32_t  v,
int32_t  placementIndex,
int32_t  type,
int32_t  partialType,
bool  reco 
)

Definition at line 883 of file HGCalCellOffset.cc.

References custom_jme_cff::area, cellArea, cellAreaPartial, cellAreaUV(), HGCalCell::extendedCell, hgcalcell_, HGCalParameters::k_ScaleToDDD2, ncell_, HGCalCell::partiaclCellsPosOffset, HGCalCell::partiaclWaferCellsOffset, HGCalCell::truncatedCell, and findQualityFiles::v.

884  {
885  if (type != 0)
886  type = 1;
887  double area(0);
888  area = cellAreaUV(u, v, placementIndex, type, reco);
889  std::pair<int, int> cell = hgcalcell_->cellType(u, v, ncell_[type], placementIndex, partialType);
890  int cellPos = cell.first;
891  int cellType = cell.second;
893  if (cellType == HGCalCell::truncatedCell || cellType == HGCalCell::extendedCell) {
894  area = reco ? cellArea[type][cellType] : HGCalParameters::k_ScaleToDDD2 * cellArea[type][cellType];
895  } else {
898  }
899  }
900  return area;
901 }
static constexpr int32_t partiaclCellsPosOffset
Definition: HGCalCell.h:70
int32_t ncell_[2]
double cellArea[2][6]
double cellAreaPartial[2][11]
static constexpr int32_t partiaclWaferCellsOffset
Definition: HGCalCell.h:39
double cellAreaUV(int32_t u, int32_t v, int32_t placementIndex, int32_t type, bool reco)
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
static constexpr double k_ScaleToDDD2

◆ cellOffsetUV2XY1() [1/2]

std::pair< double, double > HGCalCellOffset::cellOffsetUV2XY1 ( int32_t  u,
int32_t  v,
int32_t  placementIndex,
int32_t  type 
)

Definition at line 788 of file HGCalCellOffset.cc.

References HGCalCell::bottomCorner, HGCalCell::bottomLeftEdge, HGCalCell::cellPlacementExtra, HGCalCell::cornerCell, HGCalCell::extendedCell, HGCalCell::extendedMBCell, hgcalcell_, ncell_, offsetX, offsetY, HGCalCell::truncatedCell, HGCalCell::truncatedMBCell, and findQualityFiles::v.

Referenced by cellOffsetUV2XY1().

788  {
789  if (type != 0)
790  type = 1;
791  double x_off(0), y_off(0);
792  std::pair<int, int> cell = hgcalcell_->cellType(u, v, ncell_[type], placementIndex);
793  int cellPos = cell.first;
794  int cellType = cell.second;
795  if (cellType == HGCalCell::truncatedCell || cellType == HGCalCell::extendedCell) {
796  x_off = offsetX[type][cellType][cellPos - HGCalCell::bottomLeftEdge];
797  y_off = offsetY[type][cellType][cellPos - HGCalCell::bottomLeftEdge];
798  } else if ((cellType == HGCalCell::cornerCell) || (cellType == HGCalCell::truncatedMBCell) ||
799  (cellType == HGCalCell::extendedMBCell)) {
800  // The offset fo corner cells, is flipped along y-axis for 60 degree rotation of wafer
801  // and from forward to backward wafers
802  if (((placementIndex >= HGCalCell::cellPlacementExtra) && (placementIndex % 2 == 0)) ||
803  ((placementIndex < HGCalCell::cellPlacementExtra) && (placementIndex % 2 == 1))) {
804  cellPos = HGCalCell::bottomCorner + (6 + HGCalCell::bottomCorner - cellPos) % 6;
805  }
806  x_off = offsetX[type][cellType][cellPos - HGCalCell::bottomCorner];
807  y_off = offsetY[type][cellType][cellPos - HGCalCell::bottomCorner];
808  if (((placementIndex >= HGCalCell::cellPlacementExtra) && (placementIndex % 2 == 0)) ||
809  ((placementIndex < HGCalCell::cellPlacementExtra) && (placementIndex % 2 == 1))) {
810  x_off = -1 * x_off;
811  }
812  }
813  return std::make_pair(x_off, y_off);
814 }
int32_t ncell_[2]
static constexpr int32_t extendedMBCell
Definition: HGCalCell.h:33
static constexpr int32_t truncatedMBCell
Definition: HGCalCell.h:32
static constexpr int32_t bottomLeftEdge
Definition: HGCalCell.h:53
std::array< std::array< std::array< double, 6 >, 6 >, 2 > offsetY
static constexpr int32_t truncatedCell
Definition: HGCalCell.h:30
static constexpr int32_t extendedCell
Definition: HGCalCell.h:31
std::unique_ptr< HGCalCell > hgcalcell_
static constexpr int32_t bottomCorner
Definition: HGCalCell.h:59
static constexpr int32_t cornerCell
Definition: HGCalCell.h:29
std::array< std::array< std::array< double, 6 >, 6 >, 2 > offsetX
static constexpr int32_t cellPlacementExtra
Definition: HGCalCell.h:24

◆ cellOffsetUV2XY1() [2/2]

std::pair< double, double > HGCalCellOffset::cellOffsetUV2XY1 ( int32_t  u,
int32_t  v,
int32_t  placementIndex,
int32_t  type,
int32_t  partialType 
)

Definition at line 816 of file HGCalCellOffset.cc.

References HGCalCell::bottomCell, HGCalCell::bottomLeftEdge, cellOffsetUV2XY1(), HGCalCell::cellPlacementExtra, HGCalCell::extendedCell, HGCalCell::halfCell, HGCalCell::halfExtCell, HGCalCell::halfTrunCell, HGCalCell::HDPartial0920Cell, HGCalCell::HDPartial1021Cell, hgcalcell_, HGCalCell::LDPartial0007Cell, HGCalCell::LDPartial0209Cell, HGCalCell::LDPartial0714Cell, HGCalCell::LDPartial0815Cell, HGCalCell::LDPartial1415Cell, HGCalCell::LDPartial1515Cell, HGCalCell::leftCell, ncell_, HLT_IsoTrack_cff::offset, offsetPartialX, offsetPartialY, offsetX, offsetY, HGCalCell::partiaclCellsPosOffset, HGCalCell::partiaclWaferCellsOffset, HGCalCell::rightCell, HGCalCell::topCell, HGCalCell::topRightEdge, HGCalCell::truncatedCell, and findQualityFiles::v.

817  {
818  if (type != 0)
819  type = 1;
820  std::pair<double, double> offset = HGCalCellOffset::cellOffsetUV2XY1(u, v, placementIndex, type);
821  double x_off = offset.first;
822  double y_off = offset.second;
823  std::pair<int, int> cell = hgcalcell_->cellType(u, v, ncell_[type], placementIndex, partialType);
824  int cellPos = cell.first;
825  int cellType = cell.second;
827  if (cellType == HGCalCell::truncatedCell || cellType == HGCalCell::extendedCell) {
828  if (cellPos == HGCalCell::topCell) {
829  int Pos(0);
831  x_off = (placementIndex >= HGCalCell::cellPlacementExtra) ? -offsetX[type][cellType][Pos]
832  : offsetX[type][cellType][Pos];
833  y_off = offsetY[type][cellType][Pos];
834  } else if (cellPos == HGCalCell::bottomCell) {
835  int Pos(0);
836  Pos = (placementIndex) % HGCalCell::cellPlacementExtra;
837  x_off = (placementIndex >= HGCalCell::cellPlacementExtra) ? -offsetX[type][cellType][Pos]
838  : offsetX[type][cellType][Pos];
839  y_off = offsetY[type][cellType][Pos];
840  }
841  } else if ((cellType == HGCalCell::halfCell) || (cellType == HGCalCell::halfTrunCell) ||
842  (cellType == HGCalCell::halfExtCell) || (cellType == HGCalCell::LDPartial0714Cell) ||
843  (cellType == HGCalCell::LDPartial0815Cell) || (cellType == HGCalCell::HDPartial0920Cell) ||
844  (cellType == HGCalCell::HDPartial1021Cell)) {
845  int cellType1 = cellType - HGCalCell::partiaclWaferCellsOffset;
846  if (cellType == HGCalCell::halfCell) {
847 #ifdef EDM_ML_DEBUG
848  edm::LogVerbatim("HGCalGeom") << u << ":" << v << " 2";
849 #endif
850  }
851  if (cellPos == HGCalCell::leftCell) {
852  int placeIndex = placementIndex % HGCalCell::cellPlacementExtra;
853  x_off = offsetPartialX[type][cellType1][placeIndex];
854  y_off = offsetPartialY[type][cellType1][placeIndex];
855  } else if (cellPos == HGCalCell::rightCell) {
856  int placeIndex = (HGCalCell::cellPlacementExtra - placementIndex) % HGCalCell::cellPlacementExtra;
857  x_off = -offsetPartialX[type][cellType1][placeIndex];
858  y_off = offsetPartialY[type][cellType1][placeIndex];
859  }
860  x_off = placementIndex < HGCalCell::cellPlacementExtra ? x_off : -x_off;
861  } else if ((cellType == HGCalCell::LDPartial0209Cell) || (cellType == HGCalCell::LDPartial0007Cell) ||
862  (cellType == HGCalCell::LDPartial1415Cell) || (cellType == HGCalCell::LDPartial1515Cell)) {
863  int cellType1 = cellType - HGCalCell::partiaclWaferCellsOffset;
864  int placeIndex = placementIndex % HGCalCell::cellPlacementExtra;
865  x_off = offsetPartialX[type][cellType1][placeIndex];
866  y_off = offsetPartialY[type][cellType1][placeIndex];
867  x_off = placementIndex < HGCalCell::cellPlacementExtra ? x_off : -x_off;
868  }
869  }
870  return std::make_pair(x_off, y_off);
871 }
static constexpr int32_t LDPartial0007Cell
Definition: HGCalCell.h:43
Log< level::Info, true > LogVerbatim
static constexpr int32_t partiaclCellsPosOffset
Definition: HGCalCell.h:70
int32_t ncell_[2]
std::pair< double, double > cellOffsetUV2XY1(int32_t u, int32_t v, int32_t placementIndex, int32_t type)
std::array< std::array< std::array< double, 6 >, 11 >, 2 > offsetPartialX
static constexpr int32_t HDPartial1021Cell
Definition: HGCalCell.h:49
static constexpr int32_t HDPartial0920Cell
Definition: HGCalCell.h:48
static constexpr int32_t partiaclWaferCellsOffset
Definition: HGCalCell.h:39
static constexpr int32_t leftCell
Definition: HGCalCell.h:66
static constexpr int32_t LDPartial0209Cell
Definition: HGCalCell.h:42
static constexpr int32_t topCell
Definition: HGCalCell.h:68
static constexpr int32_t halfCell
Definition: HGCalCell.h:36
static constexpr int32_t LDPartial0714Cell
Definition: HGCalCell.h:41
static constexpr int32_t LDPartial1515Cell
Definition: HGCalCell.h:46
static constexpr int32_t bottomLeftEdge
Definition: HGCalCell.h:53
std::array< std::array< std::array< double, 6 >, 6 >, 2 > offsetY
static constexpr int32_t truncatedCell
Definition: HGCalCell.h:30
static constexpr int32_t halfTrunCell
Definition: HGCalCell.h:37
static constexpr int32_t extendedCell
Definition: HGCalCell.h:31
static constexpr int32_t topRightEdge
Definition: HGCalCell.h:56
std::unique_ptr< HGCalCell > hgcalcell_
static constexpr int32_t LDPartial0815Cell
Definition: HGCalCell.h:44
static constexpr int32_t LDPartial1415Cell
Definition: HGCalCell.h:45
static constexpr int32_t halfExtCell
Definition: HGCalCell.h:38
static constexpr int32_t rightCell
Definition: HGCalCell.h:67
std::array< std::array< std::array< double, 6 >, 11 >, 2 > offsetPartialY
std::array< std::array< std::array< double, 6 >, 6 >, 2 > offsetX
static constexpr int32_t bottomCell
Definition: HGCalCell.h:69
static constexpr int32_t cellPlacementExtra
Definition: HGCalCell.h:24

Member Data Documentation

◆ cellArea

double HGCalCellOffset::cellArea[2][6]
private

Definition at line 31 of file HGCalCellOffset.h.

Referenced by cellAreaUV(), and HGCalCellOffset().

◆ cellAreaPartial

double HGCalCellOffset::cellAreaPartial[2][11]
private

Definition at line 31 of file HGCalCellOffset.h.

Referenced by cellAreaUV(), and HGCalCellOffset().

◆ cellX_

double HGCalCellOffset::cellX_[2]
private

Definition at line 31 of file HGCalCellOffset.h.

Referenced by HGCalCellOffset().

◆ cellY_

double HGCalCellOffset::cellY_[2]
private

Definition at line 31 of file HGCalCellOffset.h.

Referenced by HGCalCellOffset().

◆ fullArea

double HGCalCellOffset::fullArea[2]
private

Definition at line 31 of file HGCalCellOffset.h.

◆ hgcalcell_

std::unique_ptr<HGCalCell> HGCalCellOffset::hgcalcell_
private

Definition at line 32 of file HGCalCellOffset.h.

Referenced by cellAreaUV(), cellOffsetUV2XY1(), and HGCalCellOffset().

◆ ncell_

int32_t HGCalCellOffset::ncell_[2]
private

Definition at line 30 of file HGCalCellOffset.h.

Referenced by cellAreaUV(), cellOffsetUV2XY1(), and HGCalCellOffset().

◆ offsetPartialX

std::array<std::array<std::array<double, 6>, 11>, 2> HGCalCellOffset::offsetPartialX
private

Definition at line 29 of file HGCalCellOffset.h.

Referenced by cellOffsetUV2XY1(), and HGCalCellOffset().

◆ offsetPartialY

std::array<std::array<std::array<double, 6>, 11>, 2> HGCalCellOffset::offsetPartialY
private

Definition at line 29 of file HGCalCellOffset.h.

Referenced by cellOffsetUV2XY1(), and HGCalCellOffset().

◆ offsetX

std::array<std::array<std::array<double, 6>, 6>, 2> HGCalCellOffset::offsetX
private

Definition at line 28 of file HGCalCellOffset.h.

Referenced by cellOffsetUV2XY1(), and HGCalCellOffset().

◆ offsetY

std::array<std::array<std::array<double, 6>, 6>, 2> HGCalCellOffset::offsetY
private

Definition at line 28 of file HGCalCellOffset.h.

Referenced by cellOffsetUV2XY1(), and HGCalCellOffset().

◆ sqrt3_

const double HGCalCellOffset::sqrt3_ = std::sqrt(3.0)
private

Definition at line 26 of file HGCalCellOffset.h.

Referenced by HGCalCellOffset().

◆ sqrt3By2_

const double HGCalCellOffset::sqrt3By2_ = (0.5 * sqrt3_)
private

Definition at line 27 of file HGCalCellOffset.h.

Referenced by HGCalCellOffset().