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, dqmdumpme::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  //std::cout << "h1 " << h1 << " h2 " << h2 << " H " << H << " cutarea1 " << cutArea1 << " cutarea2 " << cutArea2 << " cutarea3 " << cutArea3 << " " << cellX_[k] << " " << guardRingSizeOffset_ << std::endl;
47  double x3_1 = -(((2.0 * std::pow(h, 3)) / (3.0 * sqrt3_) - cellX_[k] * std::pow(h, 2)) / A1);
48  double y3_1 = 0;
49  double x3_2 = -(sqrt3By2_ * cellX_[k] / 3);
50  double y3_2 = cellX_[k] / 6.0;
51  double x3_3 = -(cellX_[k] * sqrt3_ / 4.0);
52  double y3_3 = cellX_[k] * 11.0 / 12.0;
53 
54  double x1 = -(3.0 * sqrt3_ * cellX_[k] / 8.0 - sqrt3_ * guardRingOffset_ / 4.0);
55  double y1 = 5.0 * cellX_[k] / 12.0 - guardRingOffset_ / 4.0;
56  double x2 = -((0.5 * cellX_[k] - 0.5 * guardRingOffset_) * sqrt3By2_);
57  double y2 = -(0.5 * cellX_[k] - 0.5 * guardRingOffset_) * 0.5;
58  double x3 = (A1 * x3_1 + A2 * x3_2 + A3 * x3_3) / cutArea3;
59  double y3 = (A1 * y3_1 + A2 * y3_2 + A3 * y3_3) / cutArea3;
60  cellArea[k][j] = totalArea - cutArea1 - cutArea2 - cutArea3;
61  double xMag1 =
62  ((5.0 * sqrt3_ * cellX_[k] / 132.0) * totalArea - (cutArea1 * x1) - (cutArea2 * x2) - (cutArea3 * x3)) /
63  (cellArea[k][j]);
64  double yMag1 =
65  ((19.0 * cellX_[k] / 132.0) * totalArea - (cutArea1 * y1) - (cutArea2 * y2) - (cutArea3 * y3)) /
66  (cellArea[k][j]);
67 
68  double xMag = 0.5 * xMag1 + sqrt3By2_ * yMag1;
69  double yMag = sqrt3By2_ * xMag1 - 0.5 * yMag1;
70 
71  std::array<double, 6> tempOffsetX = {{(sqrt3By2_ * xMag - 0.5 * yMag),
72  yMag,
73  yMag,
74  (sqrt3By2_ * xMag - 0.5 * yMag),
75  (-sqrt3By2_ * xMag - 0.5 * yMag),
76  (-sqrt3By2_ * xMag - 0.5 * yMag)}};
77  std::array<double, 6> tempOffsetY = {{(0.5 * xMag + sqrt3By2_ * yMag),
78  xMag,
79  -xMag,
80  (-0.5 * xMag - sqrt3By2_ * yMag),
81  (0.5 * xMag - sqrt3By2_ * yMag),
82  (-0.5 * xMag + sqrt3By2_ * yMag)}};
83 
84  for (int i = 0; i < 6; ++i) {
85  offsetX[k][j][i] = tempOffsetX[i];
86  offsetY[k][j][i] = tempOffsetY[i];
87  }
88  } else if (k == 1) {
89  double h = (mouseBiteCut_ - guardRingOffset_) / sqrt3By2_ - cellX_[k] / 2;
90  double H = mouseBiteCut_ + guardRingOffset_ - (1 / sqrt3By2_ * guardRingOffset_);
91  double h1 = H - (sqrt3_ / 4 * cellX_[k]) + (guardRingSizeOffset_ / (2 * sqrt3_));
92  double h2 = H - (sqrt3_ / 2 * cellX_[k]) + (guardRingSizeOffset_ / (2 * sqrt3_));
93  double totalArea = 11.0 * sqrt3_ * std::pow(cellX_[k], 2) / 8.0;
94  double cutArea1 =
95  (sqrt3By2_ * cellX_[k] * guardRingSizeOffset_) - (0.5 / sqrt3_ * std::pow(guardRingSizeOffset_, 2));
96  double cutArea2 =
97  (sqrt3_ * cellX_[k] * guardRingSizeOffset_) - (0.5 / sqrt3_ * std::pow(guardRingSizeOffset_, 2));
98  double cutArea3 =
99  sqrt3_ * std::pow(H, 2) - (1 / sqrt3By2_ * std::pow(h1, 2)) - (1 / sqrt3By2_ * std::pow(h2, 2));
100  //double cutArea3 = sqrt3_ * std::pow((mouseBiteCut_ - guardRingOffset_), 2) - sqrt3By2_ * std::pow(h, 2);
101 
102  double x2_0 = (0.375 * cellX_[k] * cellX_[k] - (0.25 * cellX_[k] * guardRingOffset_) +
103  std::pow(guardRingOffset_, 2) / 18) /
104  (sqrt3By2_ * cellX_[k] + guardRingOffset_ / (2 * sqrt3_));
105  double y2_0 = (sqrt3_ * cellX_[k] * guardRingOffset_ / 4 + std::pow(guardRingOffset_, 2) / (6 * sqrt3_)) /
106  (sqrt3By2_ * cellX_[k] + guardRingOffset_ / (2 * sqrt3_));
107  double x3_1 = -(cellX_[k] - guardRingOffset_ - 2 * (mouseBiteCut_ - guardRingOffset_) / 3) * sqrt3By2_;
108  double y3_1 = 0.5 * (cellX_[k] - guardRingOffset_ - 2 * (mouseBiteCut_ - guardRingOffset_) / 3);
109  double x3_2 = -((3 * cellX_[k] / 2 - h / 3) * sqrt3By2_ + sqrt3_ * h / 6);
110  double y3_2 = -(cellX_[k] / 4 + 4 * h / 6);
111  double A1 = sqrt3_ * std::pow((mouseBiteCut_ - guardRingOffset_), 2);
112  double A2 = sqrt3By2_ * std::pow(h, 2);
113 
114  double x1 = 0;
115  double y1 = 0.5 * cellX_[k] - 0.5 * guardRingOffset_;
116  double x2 = -(1.5 * sqrt3By2_ * cellX_[k] - x2_0 * 0.5 - y2_0 * sqrt3By2_);
117  double y2 = -(0.25 * cellX_[k] - x2_0 * sqrt3By2_ + y2_0 / 2);
118  double x3 = (A1 * x3_1 - A2 * x3_2) / (A1 - A2);
119  double y3 = (A1 * y3_1 - A2 * y3_2) / (A1 - A2);
120  cellArea[k][j] = totalArea - cutArea1 - cutArea2 - cutArea3;
121  double xMag1 = ((0.0) * totalArea - (cutArea1 * x1) - (cutArea2 * x2) - (cutArea3 * x3)) / (cellArea[k][j]);
122  double yMag1 = ((-5 * cellX_[k] / 42) * totalArea - (cutArea1 * y1) - (cutArea2 * y2) - (cutArea3 * y3)) /
123  (cellArea[k][j]);
124 
125  double xMag = -0.5 * xMag1 - sqrt3By2_ * yMag1;
126  double yMag = sqrt3By2_ * xMag1 - 0.5 * yMag1;
127 
128  std::array<double, 6> tempOffsetX = {{(sqrt3By2_ * xMag - 0.5 * yMag),
129  yMag,
130  yMag,
131  (sqrt3By2_ * xMag - 0.5 * yMag),
132  (-sqrt3By2_ * xMag - 0.5 * yMag),
133  (-sqrt3By2_ * xMag - 0.5 * yMag)}};
134  std::array<double, 6> tempOffsetY = {{(0.5 * xMag + sqrt3By2_ * yMag),
135  xMag,
136  -xMag,
137  (-0.5 * xMag - sqrt3By2_ * yMag),
138  (0.5 * xMag - sqrt3By2_ * yMag),
139  (-0.5 * xMag + sqrt3By2_ * yMag)}};
140  for (int i = 0; i < 6; ++i) {
141  offsetX[k][j][i] = tempOffsetX[i];
142  offsetY[k][j][i] = tempOffsetY[i];
143  }
144  }
145  } else if (j == HGCalCell::truncatedCell) { // Offset for truncated cells
146  double totalArea = (5.0 * sqrt3_ / 4.0) * std::pow(cellX_[k], 2); // Area of cell without any dead zone
147  double cutArea =
148  cellX_[k] * sqrt3_ * guardRingSizeOffset_; // Area of inactive region form guardring and other effects
149  cellArea[k][j] = totalArea - cutArea;
150  double offMag = (((-2.0 / 15.0) * totalArea * cellX_[k]) - ((cellX_[k] - (0.5 * guardRingOffset_)) * cutArea)) /
151  (cellArea[k][j]); // Magnitude of offset
152  // (x, y) coordinates of offset for 6 sides of wafer starting from bottom left edge in clockwise direction
153  // offset_x = -Offset_magnitude * sin(30 + 60*i) i in (0-6)
154  // offset_y = -Offset_magnitude * cos(30 + 60*i) i in (0-6)
155  std::array<double, 6> tempOffsetX = {
156  {-0.5 * offMag, -offMag, -0.5 * offMag, 0.5 * offMag, offMag, 0.5 * offMag}};
157  std::array<double, 6> tempOffsetY = {
158  {-sqrt3By2_ * offMag, 0.0, sqrt3By2_ * offMag, sqrt3By2_ * offMag, 0.0, -sqrt3By2_ * offMag}};
159  for (int i = 0; i < 6; ++i) {
160  offsetX[k][j][i] = tempOffsetX[i];
161  offsetY[k][j][i] = tempOffsetY[i];
162  }
163  } else if (j == HGCalCell::extendedCell) { //Offset for extended cells
164  double totalArea = (7.0 * sqrt3_ / 4.0) * std::pow(cellX_[k], 2); // Area of cell without any dead zone
165  double cutArea =
166  cellX_[k] * sqrt3_ * guardRingSizeOffset_; // Area of inactive region form guardring and other effects
167  cellArea[k][j] = totalArea - cutArea;
168  double offMag = // Magnitude of offset
169  (((5.0 / 42.0) * totalArea * cellX_[k]) - ((cellX_[k] - (0.5 * guardRingOffset_))) * (cutArea)) /
170  (cellArea[k][j]);
171  // (x, y) coordinates of offset for 6 sides of wafer starting from bottom left edge in clockwise direction
172  // offset_x = -Offset_magnitude * sin(30 + 60*i) i in (0-6)
173  // offset_y = -Offset_magnitude * cos(30 + 60*i) i in (0-6)
174  std::array<double, 6> tempOffsetX = {
175  {-0.5 * offMag, -offMag, -0.5 * offMag, 0.5 * offMag, offMag, 0.5 * offMag}};
176  std::array<double, 6> tempOffsetY = {
177  {-sqrt3By2_ * offMag, 0.0, sqrt3By2_ * offMag, sqrt3By2_ * offMag, 0.0, -sqrt3By2_ * offMag}};
178  for (int i = 0; i < 6; ++i) {
179  offsetX[k][j][i] = tempOffsetX[i];
180  offsetY[k][j][i] = tempOffsetY[i];
181  }
182  } else if (j == HGCalCell::truncatedMBCell) {
183  double H = mouseBiteCut_ + guardRingOffset_ - (1 / sqrt3By2_ * guardRingOffset_);
184  double h = H - (sqrt3_ / 2 * cellX_[k]) + (guardRingSizeOffset_ / (2 * sqrt3_));
185  if (h > 0) {
186  double totalArea = 5.0 * sqrt3_ * std::pow(cellX_[k], 2) / 4.0;
187 
188  double cutArea1 = (sqrt3_ * cellX_[k] * guardRingSizeOffset_);
189  double cutArea2 = std::pow(h, 2) / sqrt3By2_;
190 
191  double x1 = -(0.5 * cellX_[k] - 0.5 * guardRingOffset_) * sqrt3By2_;
192  double y1 = -(0.5 * cellX_[k] - 0.5 * guardRingOffset_) * 0.5;
193  double x2 = -((sqrt3By2_ * cellX_[k]) - (2.0 * h) / 3.0);
194  double y2 = 0.5 * cellX_[k] - (2.0 * h) / (3.0 * sqrt3_);
195  cellArea[k][j] = totalArea - cutArea1 - cutArea2;
196  //std::cout << "trunMB h " << h << " tot " << totalArea << " cutArea1 " << cutArea1 << " cutArea2 " << cutArea2 << std::endl;
197  double xMag1 =
198  ((sqrt3_ * cellX_[k] / 15.0) * totalArea - (cutArea1 * x1) - (cutArea2 * x2)) / (cellArea[k][j]);
199  double yMag1 = ((cellX_[k] / 15.0) * totalArea - (cutArea1 * y1) - (cutArea2 * y2)) / (cellArea[k][j]);
200  double xMag = -yMag1;
201  double yMag = -xMag1;
202 
203  std::array<double, 6> tempOffsetX = {{(sqrt3By2_ * xMag - 0.5 * yMag),
204  (-sqrt3By2_ * xMag - 0.5 * yMag),
205  (-sqrt3By2_ * xMag - 0.5 * yMag),
206  (sqrt3By2_ * xMag - 0.5 * yMag),
207  yMag,
208  yMag}};
209  std::array<double, 6> tempOffsetY = {{(0.5 * xMag + sqrt3By2_ * yMag),
210  (0.5 * xMag - sqrt3By2_ * yMag),
211  (-0.5 * xMag + sqrt3By2_ * yMag),
212  (-0.5 * xMag - sqrt3By2_ * yMag),
213  xMag,
214  -xMag}};
215  for (int i = 0; i < 6; ++i) {
216  offsetX[k][j][i] = tempOffsetX[i];
217  offsetY[k][j][i] = tempOffsetY[i];
218  }
219  } else {
221  std::array<double, 6> tempOffsetX = {{offsetX[k][HGCalCell::truncatedCell][0],
227  std::array<double, 6> tempOffsetY = {{offsetY[k][HGCalCell::truncatedCell][0],
233  for (int i = 0; i < 6; ++i) {
234  offsetX[k][j][i] = tempOffsetX[i];
235  offsetY[k][j][i] = tempOffsetY[i];
236  }
237  }
238  } else if (j == HGCalCell::extendedMBCell) {
239  double H = mouseBiteCut_ + guardRingOffset_ - (1 / sqrt3By2_ * guardRingOffset_);
240  double h = H - (sqrt3_ / 4 * cellX_[k]) + (guardRingSizeOffset_ / (2 * sqrt3_));
241 
242  double totalArea = 7.0 * sqrt3_ * std::pow(cellX_[k], 2) / 4.0;
243  double cutArea1 = (sqrt3_ * cellX_[k] * guardRingSizeOffset_);
244  double cutArea2 = std::pow(h, 2) / sqrt3By2_;
245 
246  double x1 = -(sqrt3By2_ * cellX_[k] - sqrt3By2_ * guardRingOffset_ / 2.0);
247  double y1 = (0.5 * cellX_[k] - 0.25 * guardRingOffset_);
248  double x2 = -(sqrt3By2_ * 1.5 * cellX_[k] - h / sqrt3_);
249  double y2 = -0.25 * cellX_[k] + h / 3.0;
250  cellArea[k][j] = totalArea - cutArea1 - cutArea2;
251  //std::cout << H << "trunMB h " << h << " tot " << totalArea << " cutArea1 " << cutArea1 << " cutArea2 " << cutArea2 << std::endl;
252  double xMag1 =
253  ((-10.0 * sqrt3_ * cellX_[k] / 168.0) * totalArea - (cutArea1 * x1) - (cutArea2 * x2)) / (cellArea[k][j]);
254  double yMag1 = ((10.0 * cellX_[k] / 168.0) * totalArea - (cutArea1 * y1) - (cutArea2 * y2)) / (cellArea[k][j]);
255 
256  double xMag = yMag1;
257  double yMag = -xMag1;
258 
259  std::array<double, 6> tempOffsetX = {{(sqrt3By2_ * xMag - 0.5 * yMag),
260  yMag,
261  yMag,
262  (sqrt3By2_ * xMag - 0.5 * yMag),
263  (-sqrt3By2_ * xMag - 0.5 * yMag),
264  (-sqrt3By2_ * xMag - 0.5 * yMag)}};
265  std::array<double, 6> tempOffsetY = {{(0.5 * xMag + sqrt3By2_ * yMag),
266  xMag,
267  -xMag,
268  (-0.5 * xMag - sqrt3By2_ * yMag),
269  (0.5 * xMag - sqrt3By2_ * yMag),
270  (-0.5 * xMag + sqrt3By2_ * yMag)}};
271  for (int i = 0; i < 6; ++i) {
272  offsetX[k][j][i] = tempOffsetX[i];
273  offsetY[k][j][i] = tempOffsetY[i];
274  }
275  }
276  }
278  ++j) { //For cells in partial wafers
279  if (j == (HGCalCell::halfCell)) {
280  double totalArea = (3.0 * sqrt3_ / 4.0) * std::pow(cellX_[k], 2);
281  double cutArea = cellX_[k] * 2.0 * guardRingOffset_ - std::pow(guardRingOffset_, 2) / sqrt3_;
282  cellAreaPartial[k][j - HGCalCell::partiaclWaferCellsOffset] = totalArea - cutArea;
283  double x1 = (-cellX_[k] * guardRingOffset_ + 2 * std::pow(guardRingOffset_, 2) / (3 * sqrt3_)) /
284  (2 * cellX_[k] - guardRingOffset_ / sqrt3_);
285  double y1 = 0;
286  double xMag = ((-2.0 * sqrt3_ * cellX_[k] / 9.0) * totalArea - (cutArea * x1)) /
288  double yMag = (0 * totalArea - (cutArea * y1)) / (cellAreaPartial[k][j - HGCalCell::partiaclWaferCellsOffset]);
289 
290  std::array<double, 6> tempOffsetX = {{(-sqrt3By2_ * xMag - 0.5 * yMag),
291  (-sqrt3By2_ * xMag + 0.5 * yMag),
292  yMag,
293  (sqrt3By2_ * xMag + 0.5 * yMag),
294  (sqrt3By2_ * xMag - 0.5 * yMag),
295  -yMag}};
296  std::array<double, 6> tempOffsetY = {{(0.5 * xMag - sqrt3By2_ * yMag),
297  (-sqrt3By2_ * yMag - 0.5 * xMag),
298  -xMag,
299  (-0.5 * xMag + sqrt3By2_ * yMag),
300  (0.5 * xMag + sqrt3By2_ * yMag),
301  xMag}};
302 
303  for (int i = 0; i < 6; ++i) {
306  }
307  } else if (j == (HGCalCell::halfTrunCell)) {
308  double totalArea = 5 * sqrt3_ * std::pow(cellX_[k], 2) / 8;
309  double cutArea1 = (sqrt3By2_ * cellX_[k] * guardRingSizeOffset_) - guardRingOffset_ * guardRingSizeOffset_;
310  double cutArea2 = (3 * cellX_[k] * guardRingOffset_) / 2 - std::pow(guardRingOffset_, 2) / (2 * sqrt3_);
311 
312  double x1 = -sqrt3_ * cellX_[k] / 4;
313  double y1 = (0.5 * cellX_[k] - 0.5 * guardRingOffset_);
314  double x2 = (-3 * cellX_[k] * guardRingOffset_ / 4 + std::pow(guardRingOffset_, 2) / (3 * sqrt3_)) /
315  (3 * cellX_[k] / 2 - guardRingOffset_ / (2 * sqrt3_));
316  double y2 = (-cellX_[k] * guardRingOffset_ / (2 * sqrt3_) + std::pow(guardRingOffset_, 2) / 18 -
317  3 * std::pow(cellX_[k], 2) / 8) /
318  (3 * cellX_[k] / 2 - guardRingOffset_ / (2 * sqrt3_));
319  cellAreaPartial[k][j - HGCalCell::partiaclWaferCellsOffset] = totalArea - cutArea1 - cutArea2;
320  double xMag1 = ((-7 * sqrt3_ * cellX_[k] / 30) * totalArea - (cutArea1 * x1) - (cutArea2 * x2)) /
322  double yMag = ((-2 * cellX_[k] / 15) * totalArea - (cutArea1 * y1) - (cutArea2 * y2)) /
324  double xMag = -xMag1;
325 
326  std::array<double, 6> tempOffsetX = {{(-sqrt3By2_ * xMag - 0.5 * yMag),
327  (-sqrt3By2_ * xMag + 0.5 * yMag),
328  yMag,
329  (sqrt3By2_ * xMag + 0.5 * yMag),
330  (sqrt3By2_ * xMag - 0.5 * yMag),
331  -yMag}};
332  std::array<double, 6> tempOffsetY = {{(0.5 * xMag - sqrt3By2_ * yMag),
333  (-sqrt3By2_ * yMag - 0.5 * xMag),
334  -xMag,
335  (-0.5 * xMag + sqrt3By2_ * yMag),
336  (0.5 * xMag + sqrt3By2_ * yMag),
337  xMag}};
338  for (int i = 0; i < 6; ++i) {
341  }
342  } else if (j == (HGCalCell::halfExtCell)) {
343  double totalArea = (7.0 * sqrt3_ / 8.0) * std::pow(cellX_[k], 2);
344  double cutArea1 = cellX_[k] * sqrt3By2_ * guardRingOffset_ - std::pow(guardRingOffset_, 2);
345  double cutArea2 = cellX_[k] * 2.0 * guardRingOffset_ - std::pow(guardRingOffset_, 2) / (2 * sqrt3_);
346  cellAreaPartial[k][j - HGCalCell::partiaclWaferCellsOffset] = totalArea - cutArea1 - cutArea2;
347 
348  double x1 = -sqrt3By2_ * cellX_[k] / 2;
349  double y1 = -(cellX_[k] - guardRingOffset_ / 2);
350  double x2 = (-cellX_[k] * guardRingOffset_ + std::pow(guardRingOffset_, 2) / (3 * sqrt3_)) /
351  (2 * cellX_[k] - guardRingOffset_ / (2 * sqrt3_));
352  double y2 = (-cellX_[k] * guardRingOffset_ / (2 * sqrt3_) + std::pow(guardRingOffset_, 2) / 18) /
353  (2 * cellX_[k] - guardRingOffset_ / (2 * sqrt3_));
354  double xMag = ((-5 * sqrt3_ * cellX_[k] / 21.0) * totalArea - (cutArea1 * x1) - (cutArea2 * x2)) /
356  double yMag = ((-5 * cellX_[k] / 42.0) * totalArea - (cutArea1 * y1) - (cutArea2 * y2)) /
358 
359  std::array<double, 6> tempOffsetX = {{(-sqrt3By2_ * xMag - 0.5 * yMag),
360  (-sqrt3By2_ * xMag + 0.5 * yMag),
361  yMag,
362  (sqrt3By2_ * xMag + 0.5 * yMag),
363  (sqrt3By2_ * xMag - 0.5 * yMag),
364  -yMag}};
365  std::array<double, 6> tempOffsetY = {{(0.5 * xMag - sqrt3By2_ * yMag),
366  (-sqrt3By2_ * yMag - 0.5 * xMag),
367  -xMag,
368  (-0.5 * xMag + sqrt3By2_ * yMag),
369  (0.5 * xMag + sqrt3By2_ * yMag),
370  xMag}};
371 
372  for (int i = 0; i < 6; ++i) {
375  }
376  } else if (j == (HGCalCell::LDPartial0714Cell)) {
377  if (k == 1) {
378  double totalArea = (9.0 * sqrt3_ / 4.0) * std::pow(cellX_[k], 2);
379  double cutArea1 =
380  (3 * cellX_[k] * sqrt3By2_ * guardRingOffset_) - (std::pow(guardRingOffset_, 2) / (2 * sqrt3_));
381  double cutArea2 = (3 * cellX_[k] * sqrt3By2_ * guardRingOffset_);
382  double cutArea3 = sqrt3_ * std::pow((mouseBiteCut_ - (guardRingOffset_ / sqrt3By2_)), 2) / 2;
383  double x1_0 = ((3.375 * cellX_[k] * cellX_[k]) - (cellX_[k] * 0.75 * guardRingOffset_) +
384  (std::pow(guardRingOffset_, 2) / 18)) /
385  ((3 * cellX_[k] * sqrt3By2_) - (guardRingOffset_ / (2 * sqrt3_)));
386  double y1_0 =
387  ((3 * cellX_[k] * sqrt3By2_ * guardRingOffset_ / 2) - (std::pow(guardRingOffset_, 2) / (6 * sqrt3_))) /
388  ((3 * cellX_[k] * sqrt3By2_) - (guardRingOffset_ / (2 * sqrt3_)));
389 
390  double x2_0 = (3 * sqrt3By2_ * cellX_[k] / 2);
391  double y2_0 = guardRingOffset_ / 2;
392 
393  double x1 = (cellX_[k] / 2 - guardRingOffset_) * sqrt3By2_ + x1_0 * 0.5 + y1_0 * sqrt3By2_;
394  double y1 = cellX_[k] + (cellX_[k] / 2 - guardRingOffset_) * 0.5 - x1_0 * sqrt3By2_ + y1_0 * 0.5;
395 
396  double x2 = x2_0 - sqrt3By2_ * cellX_[k];
397  double y2 = -(cellX_[k] - y2_0);
398 
399  double x3 = sqrt3_ * cellX_[k] - mouseBiteCut_ + (mouseBiteCut_ - (guardRingOffset_ / sqrt3By2_)) / 3;
400  double y3 = -(cellX_[k] - sqrt3_ * (mouseBiteCut_ - (guardRingOffset_ / sqrt3By2_)) / 3 - guardRingOffset_);
401 
402  cellAreaPartial[k][j - HGCalCell::partiaclWaferCellsOffset] = totalArea - cutArea1 - cutArea2 - cutArea3;
403  double xMag = ((sqrt3_ * cellX_[k] / 8) * totalArea - (cutArea1 * x1) - (cutArea2 * x2) - (cutArea3 * x3)) /
405  double yMag = ((-1 * cellX_[k] / 8) * totalArea - (cutArea1 * y1) - (cutArea2 * y2) - (cutArea3 * y3)) /
407 
408  std::array<double, 6> tempOffsetX = {{(-sqrt3By2_ * xMag - 0.5 * yMag),
409  (-sqrt3By2_ * xMag + 0.5 * yMag),
410  yMag,
411  (sqrt3By2_ * xMag + 0.5 * yMag),
412  (sqrt3By2_ * xMag - 0.5 * yMag),
413  -yMag}};
414  std::array<double, 6> tempOffsetY = {{(0.5 * xMag - sqrt3By2_ * yMag),
415  (-sqrt3By2_ * yMag - 0.5 * xMag),
416  -xMag,
417  (-0.5 * xMag + sqrt3By2_ * yMag),
418  (0.5 * xMag + sqrt3By2_ * yMag),
419  xMag}};
420 
421  for (int i = 0; i < 6; ++i) {
424  }
425  } else {
427  for (int i = 0; i < 6; ++i) {
430  }
431  }
432  } else if (j == (HGCalCell::LDPartial0209Cell)) {
433  if (k == 1) {
434  double totalArea = (23.0 * sqrt3_ / 8.0) * std::pow(cellX_[k], 2);
435  double cutArea1 =
436  (5 * cellX_[k] * sqrt3By2_ * guardRingOffset_) - (std::pow(guardRingOffset_, 2) / (2 * sqrt3_));
437  double cutArea2 = (4 * cellX_[k] * guardRingOffset_);
438  double cutArea3 = std::pow(mouseBiteCut_, 2) / sqrt3_;
439 
440  double x1_0 = (9.375 * cellX_[k] * cellX_[k] - (cellX_[k] * 1.25 * guardRingOffset_) +
441  (std::pow(guardRingOffset_, 2) / 18)) /
442  ((5 * cellX_[k] * sqrt3By2_) - (guardRingOffset_ / (2 * sqrt3_)));
443  double y1_0 =
444  ((5 * cellX_[k] * sqrt3By2_ * guardRingOffset_ / 2) - (std::pow(guardRingOffset_, 2) / (6 * sqrt3_))) /
445  ((5 * cellX_[k] * sqrt3By2_) - (guardRingOffset_ / (2 * sqrt3_)));
446 
447  double x1 = (1.5 * cellX_[k]) * sqrt3By2_ - x1_0 * 0.5 - y1_0 * sqrt3By2_;
448  double y1 = -0.25 * cellX_[k] + x1_0 * sqrt3By2_ - y1_0 * 0.5;
449  double x2 = -(sqrt3By2_ * cellX_[k] - 0.5 * guardRingOffset_);
450  double y2 = 1.5 * cellX_[k];
451  double x3 = -(sqrt3By2_ * cellX_[k] - mouseBiteCut_ / 3);
452  double y3 = 3.5 * cellX_[k] - (5 * mouseBiteCut_) / 3 * sqrt3_;
453 
454  cellAreaPartial[k][j - HGCalCell::partiaclWaferCellsOffset] = totalArea - cutArea1 - cutArea2 - cutArea3;
455  double xMag =
456  ((-9 * cellX_[k] / (sqrt3_ * 92)) * totalArea - (cutArea1 * x1) - (cutArea2 * x2) - (cutArea3 * x3)) /
458  double yMag =
459  ((199 * cellX_[k] / (sqrt3_ * 276)) * totalArea - (cutArea1 * y1) - (cutArea2 * y2) - (cutArea3 * y3)) /
461 
462  std::array<double, 6> tempOffsetX = {{(-sqrt3By2_ * xMag - 0.5 * yMag),
463  (-sqrt3By2_ * xMag + 0.5 * yMag),
464  yMag,
465  (sqrt3By2_ * xMag + 0.5 * yMag),
466  (sqrt3By2_ * xMag - 0.5 * yMag),
467  -yMag}};
468  std::array<double, 6> tempOffsetY = {{(0.5 * xMag - sqrt3By2_ * yMag),
469  (-sqrt3By2_ * yMag - 0.5 * xMag),
470  -xMag,
471  (-0.5 * xMag + sqrt3By2_ * yMag),
472  (0.5 * xMag + sqrt3By2_ * yMag),
473  xMag}};
474  for (int i = 0; i < 6; ++i) {
477  }
478  } else {
480  for (int i = 0; i < 6; ++i) {
483  }
484  }
485  } else if (j == (HGCalCell::LDPartial0007Cell)) {
486  if (k == 1) {
487  double totalArea = (5.0 * sqrt3_ / 4.0) * std::pow(cellX_[k], 2);
488  double cutArea1 = (cellX_[k] * guardRingOffset_);
489  double cutArea2 = (sqrt3_ * cellX_[k] * guardRingOffset_);
490  double h = cellX_[k] - (sqrt3By2_ * cellX_[k] - mouseBiteCut_) / sqrt3By2_;
491  double cutArea3 = sqrt3_ * std::pow(h, 2) / 2;
492 
493  double x1 = cellX_[k] * sqrt3By2_ - guardRingOffset_ / 2;
494  double y1 = 0;
495  double x2 = 0;
496  double y2 = 0.5 * cellX_[k] - guardRingOffset_ / 2;
497  double x3 = sqrt3By2_ * cellX_[k] - guardRingOffset_ - h / sqrt3_;
498  double y3 = 0.5 * cellX_[k] - guardRingOffset_ - h / 3;
499 
500  cellAreaPartial[k][j - HGCalCell::partiaclWaferCellsOffset] = totalArea - cutArea1 - cutArea2 - cutArea3;
501  double xMag = ((0.0) * totalArea - (cutArea1 * x1) - (cutArea2 * x2) - (cutArea3 * x3)) /
503  double yMag = ((-2 * cellX_[k] / 15) * totalArea - (cutArea1 * y1) - (cutArea2 * y2) - (cutArea3 * y3)) /
505 
506  std::array<double, 6> tempOffsetX = {{(-sqrt3By2_ * xMag - 0.5 * yMag),
507  (-sqrt3By2_ * xMag + 0.5 * yMag),
508  yMag,
509  (sqrt3By2_ * xMag + 0.5 * yMag),
510  (sqrt3By2_ * xMag - 0.5 * yMag),
511  -yMag}};
512  std::array<double, 6> tempOffsetY = {{(0.5 * xMag - sqrt3By2_ * yMag),
513  (-sqrt3By2_ * yMag - 0.5 * xMag),
514  -xMag,
515  (-0.5 * xMag + sqrt3By2_ * yMag),
516  (0.5 * xMag + sqrt3By2_ * yMag),
517  xMag}};
518  for (int i = 0; i < 6; ++i) {
521  }
522  } else {
524  for (int i = 0; i < 6; ++i) {
527  }
528  }
529  } else if (j == (HGCalCell::LDPartial0815Cell)) {
530  if (k == 1) {
531  double totalArea = sqrt3_ * std::pow(cellX_[k], 2);
532  double cutArea1 = (sqrt3_ * cellX_[k] * guardRingOffset_);
533  double cutArea2 = (sqrt3_ * cellX_[k] * guardRingOffset_) - std::pow(guardRingOffset_, 2) / (2 * sqrt3_);
534  double cutArea3 = sqrt3_ * std::pow((mouseBiteCut_ - guardRingOffset_ / sqrt3By2_), 2) / 2;
535 
536  double x2_0 = (1.5 * cellX_[k] * cellX_[k] - (0.5 * cellX_[k] * guardRingOffset_) +
537  std::pow(guardRingOffset_, 2) / 18) /
538  (sqrt3_ * cellX_[k] - guardRingOffset_ / (2 * sqrt3_));
539  double y2_0 = (sqrt3By2_ * cellX_[k] * guardRingOffset_ - std::pow(guardRingOffset_, 2) / (sqrt3_ * 3)) /
540  (sqrt3_ * cellX_[k] - guardRingOffset_ / (2 * sqrt3_));
541  double x1 = 0;
542  double y1 = 0.5 * cellX_[k] - guardRingOffset_ / 2;
543  double x2 = x2_0 * 0.5 - y2_0 * sqrt3By2_;
544  double y2 = -(cellX_[k] - (x2_0 * sqrt3By2_ + y2_0 * 0.5));
545  double x3 = sqrt3By2_ * cellX_[k] - mouseBiteCut_ + (mouseBiteCut_ - guardRingOffset_ / sqrt3By2_) / 3;
546  double y3 = cellX_[k] - (mouseBiteCut_ - guardRingOffset_ / sqrt3By2_) / sqrt3_ - guardRingOffset_;
547 
548  cellAreaPartial[k][j - HGCalCell::partiaclWaferCellsOffset] = totalArea - cutArea1 - cutArea2 - cutArea3;
549  double xMag = -((-sqrt3_ * cellX_[k] / 9) * totalArea - (cutArea1 * x1) - (cutArea2 * x2) - (cutArea3 * x3)) /
551  double yMag = ((-cellX_[k] / 15) * totalArea - (cutArea1 * y1) - (cutArea2 * y2) - (cutArea3 * y3)) /
553 
554  std::array<double, 6> tempOffsetX = {{(-sqrt3By2_ * xMag - 0.5 * yMag),
555  (-sqrt3By2_ * xMag + 0.5 * yMag),
556  yMag,
557  (sqrt3By2_ * xMag + 0.5 * yMag),
558  (sqrt3By2_ * xMag - 0.5 * yMag),
559  -yMag}};
560  std::array<double, 6> tempOffsetY = {{(0.5 * xMag - sqrt3By2_ * yMag),
561  (-sqrt3By2_ * yMag - 0.5 * xMag),
562  -xMag,
563  (-0.5 * xMag + sqrt3By2_ * yMag),
564  (0.5 * xMag + sqrt3By2_ * yMag),
565  xMag}};
566  for (int i = 0; i < 6; ++i) {
569  }
570  } else {
571  for (int i = 0; i < 6; ++i) {
575  }
576  }
577  } else if (j == (HGCalCell::LDPartial1415Cell)) {
578  if (k == 1) {
579  double totalArea = 7 * sqrt3_ * std::pow(cellX_[k], 2) / 4;
580  double cutArea1 = (3 * cellX_[k] * guardRingOffset_);
581  double cutArea2 = (2 * sqrt3_ * cellX_[k] * guardRingOffset_) - std::pow(guardRingOffset_, 2) * sqrt3By2_;
582  double cutArea3 = std::pow((mouseBiteCut_ - guardRingOffset_), 2) / sqrt3_;
583 
584  double x2_0 = (6 * cellX_[k] * cellX_[k] - std::pow(guardRingOffset_, 2)) /
585  (2 * sqrt3_ * cellX_[k] - guardRingOffset_ * sqrt3By2_);
586  double y2_0 = (sqrt3_ * cellX_[k] * guardRingOffset_ - std::pow(guardRingOffset_, 2) / sqrt3_) /
587  (2 * sqrt3_ * cellX_[k] - guardRingOffset_ * sqrt3By2_);
588  double x1 = -sqrt3By2_ * cellX_[k] + guardRingOffset_ / 2;
589  double y1 = -cellX_[k];
590  double x2 = sqrt3By2_ * cellX_[k] - x2_0 * 0.5 - y2_0 * sqrt3By2_;
591  double y2 = 0.5 * cellX_[k] - x2_0 * sqrt3By2_ + y2_0 * 0.5;
592  double h = (mouseBiteCut_ - guardRingOffset_) / sqrt3By2_;
593  double x3 = -(cellX_[k] - h / 3 - guardRingOffset_) * sqrt3By2_;
594  double y3 = 5 * h / 6 - 5 * cellX_[k] / 2;
595 
596  cellAreaPartial[k][j - HGCalCell::partiaclWaferCellsOffset] = totalArea - cutArea1 - cutArea2 - cutArea3;
597  double xMag =
598  ((-2 * cellX_[k] / (7 * sqrt3_)) * totalArea - (cutArea1 * x1) - (cutArea2 * x2) - (cutArea3 * x3)) /
600  double yMag = ((-cellX_[k] / 3) * totalArea - (cutArea1 * y1) - (cutArea2 * y2) - (cutArea3 * y3)) /
602 
603  std::array<double, 6> tempOffsetX = {{(-sqrt3By2_ * xMag - 0.5 * yMag),
604  (-sqrt3By2_ * xMag + 0.5 * yMag),
605  yMag,
606  (sqrt3By2_ * xMag + 0.5 * yMag),
607  (sqrt3By2_ * xMag - 0.5 * yMag),
608  -yMag}};
609  std::array<double, 6> tempOffsetY = {{(0.5 * xMag - sqrt3By2_ * yMag),
610  (-sqrt3By2_ * yMag - 0.5 * xMag),
611  -xMag,
612  (-0.5 * xMag + sqrt3By2_ * yMag),
613  (0.5 * xMag + sqrt3By2_ * yMag),
614  xMag}};
615  for (int i = 0; i < 6; ++i) {
618  }
619  } else {
620  for (int i = 0; i < 6; ++i) {
624  }
625  }
626  } else if (j == (HGCalCell::LDPartial1515Cell)) {
627  if (k == 1) {
628  double totalArea = 7 * sqrt3_ * std::pow(cellX_[k], 2) / 8;
629  double cutArea1 = (2 * cellX_[k] * guardRingOffset_);
630  double cutArea2 = (sqrt3By2_ * cellX_[k] * guardRingOffset_);
631  double cutArea3 = cellX_[k] * (mouseBiteCut_ - guardRingOffset_) - (sqrt3_ * cellX_[k] * cellX_[k] / 8);
632 
633  double x1 = -guardRingOffset_ / 2;
634  double y1 = 0;
635  double x2 = -(sqrt3By2_ * cellX_[k] / 2);
636  double y2 = -(cellX_[k] - 0.5 * guardRingOffset_);
637  double x3 = (cellX_[k] * cellX_[k] / 8 - sqrt3_ * cellX_[k] * (mouseBiteCut_ - guardRingOffset_) / 4) /
638  ((mouseBiteCut_ - guardRingOffset_) - sqrt3_ * cellX_[k] / 8);
639  double y3 =
640  (std::pow((mouseBiteCut_ - guardRingOffset_), 2) / sqrt3_ -
641  (1.25 * cellX_[k] * (mouseBiteCut_ - guardRingOffset_)) + 7 * sqrt3_ * cellX_[k] * cellX_[k] / 48) /
642  ((mouseBiteCut_ - guardRingOffset_) - sqrt3_ * cellX_[k] / 8);
643 
644  cellAreaPartial[k][j - HGCalCell::partiaclWaferCellsOffset] = totalArea - cutArea1 - cutArea2 - cutArea3;
645  double xMag = (-(cellX_[k] / (sqrt3_)) * totalArea - (cutArea1 * x1) - (cutArea2 * x2) - (cutArea3 * x3)) /
647  double yMag = ((-5 * cellX_[k] / 42) * totalArea - (cutArea1 * y1) - (cutArea2 * y2) - (cutArea3 * y3)) /
649 
650  std::array<double, 6> tempOffsetX = {{(-sqrt3By2_ * xMag - 0.5 * yMag),
651  (-sqrt3By2_ * xMag + 0.5 * yMag),
652  yMag,
653  (sqrt3By2_ * xMag + 0.5 * yMag),
654  (sqrt3By2_ * xMag - 0.5 * yMag),
655  -yMag}};
656  std::array<double, 6> tempOffsetY = {{(0.5 * xMag - sqrt3By2_ * yMag),
657  (-sqrt3By2_ * yMag - 0.5 * xMag),
658  -xMag,
659  (-0.5 * xMag + sqrt3By2_ * yMag),
660  (0.5 * xMag + sqrt3By2_ * yMag),
661  xMag}};
662  for (int i = 0; i < 6; ++i) {
665  }
666  } else {
667  for (int i = 0; i < 6; ++i) {
671  }
672  }
673  } else if (j == (HGCalCell::HDPartial0920Cell)) {
674  if (k == 0) {
675  double totalArea = 37 * sqrt3_ * std::pow(cellX_[k], 2) / 24;
676  double cutArea1 = (4 * cellX_[k] * guardRingOffset_) / sqrt3_;
677  double cutArea2 =
678  (7 * cellX_[k] * guardRingOffset_) / (2 * sqrt3_) - std::pow(guardRingOffset_, 2) / (2 * sqrt3_);
679 
680  double x1 = cellX_[k] / (2 * sqrt3_);
681  double y1 = -(0.5 * cellX_[k] - 0.5 * guardRingOffset_);
682  double x2_0 = ((2.041 * cellX_[k] * cellX_[k]) - (cellX_[k] * 0.583 * guardRingOffset_) +
683  (std::pow(guardRingOffset_, 2) / 18)) /
684  ((7 * cellX_[k] / (2 * sqrt3_)) - (guardRingOffset_ / (2 * sqrt3_)));
685  double y2_0 =
686  ((7 * cellX_[k] * guardRingOffset_ / (4 * sqrt3_)) - std::pow(guardRingOffset_, 2) / (6 * sqrt3_)) /
687  ((7 * cellX_[k] / (2 * sqrt3_)) - (guardRingOffset_ / (2 * sqrt3_)));
688 
689  double x2 = (0.5 * x2_0) - (sqrt3By2_ * y2_0) + (cellX_[k] * 0.5 * sqrt3By2_);
690  double y2 = -(0.5 * y2_0) - (sqrt3By2_ * x2_0) + (cellX_[k] * 1.25);
691  cellAreaPartial[k][j - HGCalCell::partiaclWaferCellsOffset] = totalArea - cutArea1 - cutArea2;
692  double xMag = ((25 * sqrt3_ * cellX_[k] / 148) * totalArea - (cutArea1 * x1) - (cutArea2 * x2)) /
694  double yMag = ((73 * cellX_[k] / 444) * totalArea - (cutArea1 * y1) - (cutArea2 * y2)) /
696 
697  std::array<double, 6> tempOffsetX = {{(-sqrt3By2_ * xMag - 0.5 * yMag),
698  (-sqrt3By2_ * xMag + 0.5 * yMag),
699  yMag,
700  (sqrt3By2_ * xMag + 0.5 * yMag),
701  (sqrt3By2_ * xMag - 0.5 * yMag),
702  -yMag}};
703  std::array<double, 6> tempOffsetY = {{(0.5 * xMag - sqrt3By2_ * yMag),
704  (-sqrt3By2_ * yMag - 0.5 * xMag),
705  -xMag,
706  (-0.5 * xMag + sqrt3By2_ * yMag),
707  (0.5 * xMag + sqrt3By2_ * yMag),
708  xMag}};
709  for (int i = 0; i < 6; ++i) {
712  }
713  } else {
714  for (int i = 0; i < 6; ++i) {
718  }
719  }
720  } else if (j == (HGCalCell::HDPartial1021Cell)) {
721  if (k == 0) {
722  double totalArea = 11 * sqrt3_ * std::pow(cellX_[k], 2) / 6;
723  double cutArea1 = (5 * cellX_[k] * guardRingOffset_) / (2 * sqrt3_);
724  double cutArea2 =
725  (5 * cellX_[k] * guardRingOffset_) / (2 * sqrt3_) - std::pow(guardRingOffset_, 2) / (2 * sqrt3_);
726 
727  double x1 = -cellX_[k] / (4 * sqrt3_);
728  double y1 = cellX_[k] - 0.5 * guardRingOffset_;
729  double x2_0 = ((1.041 * cellX_[k] * cellX_[k]) - (cellX_[k] * 0.416 * guardRingOffset_) +
730  (std::pow(guardRingOffset_, 2) / 18.0)) /
731  ((5.0 * cellX_[k] / (2.0 * sqrt3_)) - (guardRingOffset_ / (2.0 * sqrt3_)));
732  double y2_0 =
733  ((5.0 * cellX_[k] * guardRingOffset_ / (4.0 * sqrt3_)) - std::pow(guardRingOffset_, 2) / (6 * sqrt3_)) /
734  ((5.0 * cellX_[k] / (2.0 * sqrt3_)) - (guardRingOffset_ / (2.0 * sqrt3_)));
735 
736  double x2 = -(0.5 * x2_0) + (sqrt3By2_ * y2_0) + (cellX_[k] * 1.5 * sqrt3By2_);
737  double y2 = -(0.5 * y2_0) + (sqrt3By2_ * x2_0) - cellX_[k];
738  cellAreaPartial[k][j - HGCalCell::partiaclWaferCellsOffset] = totalArea - cutArea1 - cutArea2;
739  double xMag = ((47.0 * cellX_[k] / (528.0 * sqrt3_)) * totalArea - (cutArea1 * x1) - (cutArea2 * x2)) /
741  double yMag = ((47.0 * cellX_[k] / 528.0) * totalArea - (cutArea1 * y1) - (cutArea2 * y2)) /
743 
744  std::array<double, 6> tempOffsetX = {{(-sqrt3By2_ * xMag - 0.5 * yMag),
745  (-sqrt3By2_ * xMag + 0.5 * yMag),
746  yMag,
747  (sqrt3By2_ * xMag + 0.5 * yMag),
748  (sqrt3By2_ * xMag - 0.5 * yMag),
749  -yMag}};
750  std::array<double, 6> tempOffsetY = {{(0.5 * xMag - sqrt3By2_ * yMag),
751  (-sqrt3By2_ * yMag - 0.5 * xMag),
752  -xMag,
753  (-0.5 * xMag + sqrt3By2_ * yMag),
754  (0.5 * xMag + sqrt3By2_ * yMag),
755  xMag}};
756  for (int i = 0; i < 6; ++i) {
759  }
760  } else {
761  for (int i = 0; i < 6; ++i) {
765  }
766  }
767  }
768  }
769  }
770 
771 #ifdef EDM_ML_DEBUG
772  edm::LogVerbatim("HGCalGeom") << "HGCalCellOffset initialized with waferSize " << waferSize << " number of cells "
773  << nFine << ":" << nCoarse << " Guardring offset " << guardRingOffset_ << " Mousebite "
774  << mouseBiteCut_;
775 #endif
776 }
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 861 of file HGCalCellOffset.cc.

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

Referenced by cellAreaUV().

861  {
862  if (type != 0)
863  type = 1;
864  double area(0);
865  std::pair<int, int> cell = hgcalcell_->cellType(u, v, ncell_[type], placementIndex);
866  int cellType = cell.second;
867  area = reco ? cellArea[type][cellType] : HGCalParameters::k_ScaleToDDD2 * cellArea[type][cellType];
868  return area;
869 }
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 871 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.

872  {
873  if (type != 0)
874  type = 1;
875  double area(0);
876  area = cellAreaUV(u, v, placementIndex, type, reco);
877  std::pair<int, int> cell = hgcalcell_->cellType(u, v, ncell_[type], placementIndex, partialType);
878  int cellPos = cell.first;
879  int cellType = cell.second;
881  if (cellType == HGCalCell::truncatedCell || cellType == HGCalCell::extendedCell) {
882  area = reco ? cellArea[type][cellType] : HGCalParameters::k_ScaleToDDD2 * cellArea[type][cellType];
883  } else {
886  }
887  }
888  return area;
889 }
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 778 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().

778  {
779  if (type != 0)
780  type = 1;
781  double x_off(0), y_off(0);
782  std::pair<int, int> cell = hgcalcell_->cellType(u, v, ncell_[type], placementIndex);
783  int cellPos = cell.first;
784  int cellType = cell.second;
785  if (cellType == HGCalCell::truncatedCell || cellType == HGCalCell::extendedCell) {
786  x_off = offsetX[type][cellType][cellPos - HGCalCell::bottomLeftEdge];
787  y_off = offsetY[type][cellType][cellPos - HGCalCell::bottomLeftEdge];
788  } else if ((cellType == HGCalCell::cornerCell) || (cellType == HGCalCell::truncatedMBCell) ||
789  (cellType == HGCalCell::extendedMBCell)) {
790  // The offset fo corner cells, is flipped along y-axis for 60 degree rotation of wafer
791  // and from forward to backward wafers
792  if (((placementIndex >= HGCalCell::cellPlacementExtra) && (placementIndex % 2 == 0)) ||
793  ((placementIndex < HGCalCell::cellPlacementExtra) && (placementIndex % 2 == 1))) {
794  cellPos = HGCalCell::bottomCorner + (6 + HGCalCell::bottomCorner - cellPos) % 6;
795  }
796  x_off = offsetX[type][cellType][cellPos - HGCalCell::bottomCorner];
797  y_off = offsetY[type][cellType][cellPos - HGCalCell::bottomCorner];
798  if (((placementIndex >= HGCalCell::cellPlacementExtra) && (placementIndex % 2 == 0)) ||
799  ((placementIndex < HGCalCell::cellPlacementExtra) && (placementIndex % 2 == 1))) {
800  x_off = -1 * x_off;
801  }
802  }
803  return std::make_pair(x_off, y_off);
804 }
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 806 of file HGCalCellOffset.cc.

References HGCalCell::bottomCell, HGCalCell::bottomLeftEdge, cellOffsetUV2XY1(), HGCalCell::cellPlacementExtra, gather_cfg::cout, 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_, hltrates_dqm_sourceclient-live_cfg::offset, offsetPartialX, offsetPartialY, offsetX, offsetY, HGCalCell::partiaclCellsPosOffset, HGCalCell::partiaclWaferCellsOffset, HGCalCell::rightCell, HGCalCell::topCell, HGCalCell::topRightEdge, HGCalCell::truncatedCell, and findQualityFiles::v.

807  {
808  if (type != 0)
809  type = 1;
810  std::pair<double, double> offset = HGCalCellOffset::cellOffsetUV2XY1(u, v, placementIndex, type);
811  double x_off = offset.first;
812  double y_off = offset.second;
813  std::pair<int, int> cell = hgcalcell_->cellType(u, v, ncell_[type], placementIndex, partialType);
814  int cellPos = cell.first;
815  int cellType = cell.second;
817  if (cellType == HGCalCell::truncatedCell || cellType == HGCalCell::extendedCell) {
818  if (cellPos == HGCalCell::topCell) {
819  int Pos(0);
821  x_off = (placementIndex >= HGCalCell::cellPlacementExtra) ? -offsetX[type][cellType][Pos]
822  : offsetX[type][cellType][Pos];
823  y_off = offsetY[type][cellType][Pos];
824  } else if (cellPos == HGCalCell::bottomCell) {
825  int Pos(0);
826  Pos = (placementIndex) % HGCalCell::cellPlacementExtra;
827  x_off = (placementIndex >= HGCalCell::cellPlacementExtra) ? -offsetX[type][cellType][Pos]
828  : offsetX[type][cellType][Pos];
829  y_off = offsetY[type][cellType][Pos];
830  }
831  } else if ((cellType == HGCalCell::halfCell) || (cellType == HGCalCell::halfTrunCell) ||
832  (cellType == HGCalCell::halfExtCell) || (cellType == HGCalCell::LDPartial0714Cell) ||
833  (cellType == HGCalCell::LDPartial0815Cell) || (cellType == HGCalCell::HDPartial0920Cell) ||
834  (cellType == HGCalCell::HDPartial1021Cell)) {
835  int cellType1 = cellType - HGCalCell::partiaclWaferCellsOffset;
836  if (cellType == HGCalCell::halfCell) {
837  std::cout << u << ":" << v << " 2" << std::endl;
838  }
839  if (cellPos == HGCalCell::leftCell) {
840  int placeIndex = placementIndex % HGCalCell::cellPlacementExtra;
841  x_off = offsetPartialX[type][cellType1][placeIndex];
842  y_off = offsetPartialY[type][cellType1][placeIndex];
843  } else if (cellPos == HGCalCell::rightCell) {
844  int placeIndex = (HGCalCell::cellPlacementExtra - placementIndex) % HGCalCell::cellPlacementExtra;
845  x_off = -offsetPartialX[type][cellType1][placeIndex];
846  y_off = offsetPartialY[type][cellType1][placeIndex];
847  }
848  x_off = placementIndex < HGCalCell::cellPlacementExtra ? x_off : -x_off;
849  } else if ((cellType == HGCalCell::LDPartial0209Cell) || (cellType == HGCalCell::LDPartial0007Cell) ||
850  (cellType == HGCalCell::LDPartial1415Cell) || (cellType == HGCalCell::LDPartial1515Cell)) {
851  int cellType1 = cellType - HGCalCell::partiaclWaferCellsOffset;
852  int placeIndex = placementIndex % HGCalCell::cellPlacementExtra;
853  x_off = offsetPartialX[type][cellType1][placeIndex];
854  y_off = offsetPartialY[type][cellType1][placeIndex];
855  x_off = placementIndex < HGCalCell::cellPlacementExtra ? x_off : -x_off;
856  }
857  }
858  return std::make_pair(x_off, y_off);
859 }
static constexpr int32_t LDPartial0007Cell
Definition: HGCalCell.h:43
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().